]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
ALIROOT-5559 Fix for memory jump in CPass0 (LHC10 processing). A custom streamer...
authorhristov <Peter.Hristov@cern.ch>
Wed, 13 Aug 2014 15:53:45 +0000 (17:53 +0200)
committerhristov <Peter.Hristov@cern.ch>
Wed, 13 Aug 2014 15:57:29 +0000 (17:57 +0200)
TPC/Calib/AliTPCcalibAlign.cxx
TPC/Calib/AliTPCcalibAlign.h
TPC/TPCcalibLinkDef.h

index 1bc2fff434e02643f601a620cbbf6ae371297e31..bbb60b4def3948800c36d871cb02c8d387d55405 100644 (file)
 #include "Riostream.h"
 #include "TRandom.h"
 #include <sstream>
+
+#include "AliSysInfo.h"
 using namespace std;
 
 AliTPCcalibAlign* AliTPCcalibAlign::fgInstance = 0;
@@ -3816,3 +3818,151 @@ void AliTPCcalibAlign::MakeReportDyPhi(TFile */*output*/){
 }
 
 
+//______________________________________________________________________________
+void AliTPCcalibAlign::Streamer(TBuffer &R__b)
+{
+  // Stream an object of class AliTPCcalibAlign.
+  Bool_t isDebug=AliLog::GetDebugLevel("","AliTPCcalibAlign")>0;
+  if (isDebug) AliSysInfo::SetVerbose(kTRUE);
+
+  if (R__b.IsReading()) {
+    //    R__b.ReadClassBuffer(AliTPCcalibAlign::Class(),this);
+    UInt_t R__s, R__c;
+    Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
+
+    AliTPCcalibBase::Streamer(R__b); // Stream the base class
+
+    // Read the "big data members" from the same Root directory
+    if (gDirectory){
+      for (Int_t i=0; i<2; ++i){
+       TString hisName=TString::Format("AliTPCcalibAlign.%s.fClusterDelta_%d",GetName(),i);
+       if (gDirectory->Get(hisName.Data())){
+         fClusterDelta[i] = dynamic_cast<THn*>(gDirectory->Get((hisName).Data()));
+       }
+      }
+    
+      for (Int_t i=0; i<4; ++i){
+       TString hisName=TString::Format("AliTPCcalibAlign.%s.fTrackletDelta_%d",GetName(),i);
+       if (gDirectory->Get(hisName.Data())){
+         fTrackletDelta[i] = dynamic_cast<THnSparse*>(gDirectory->Get((hisName).Data()));
+       }
+      }
+    }
+
+    // If the "big data members"were not in the file, try to read them from the object
+    // This is suppose to restore the backward compatibility
+    for (Int_t i=0; i<2; ++i) {
+      if (!fClusterDelta[i]) R__b >> fClusterDelta[i];
+    }
+
+    for (Int_t i=0; i<4; ++i) {
+      if (!fTrackletDelta[i]) R__b >> fTrackletDelta[i];
+    }
+
+    // Read all the other data members. This is error prone,
+    // please make sure this part is updated if you change the data members
+    fDphiHistArray.Streamer(R__b);
+    fDthetaHistArray.Streamer(R__b);
+    fDyHistArray.Streamer(R__b);
+    fDzHistArray.Streamer(R__b);
+
+    fDyPhiHistArray.Streamer(R__b);
+    fDzThetaHistArray.Streamer(R__b);
+    
+    fDphiZHistArray.Streamer(R__b);
+    fDthetaZHistArray.Streamer(R__b);
+    fDyZHistArray.Streamer(R__b);
+    fDzZHistArray.Streamer(R__b);
+      
+    fFitterArray12.Streamer(R__b);
+    fFitterArray9.Streamer(R__b);
+    fFitterArray6.Streamer(R__b);
+    
+    fMatrixArray12.Streamer(R__b);
+    fMatrixArray9.Streamer(R__b);
+    fMatrixArray6.Streamer(R__b);
+
+    fCombinedMatrixArray6.Streamer(R__b);
+
+    R__b.ReadFastArray(fPoints,72*72);
+
+    R__b >> fNoField;
+    R__b >> fXIO;
+    R__b >> fXmiddle;
+    R__b >> fXquadrant;
+    
+    fArraySectorIntParam.Streamer(R__b);
+    fArraySectorIntCovar.Streamer(R__b);
+
+    R__b >> fSectorParamA;
+    R__b >> fSectorCovarA;
+    R__b >> fSectorParamC;
+    R__b >> fSectorCovarC;
+    
+    R__b >> fUseInnerOuter;
+    R__b >> fgkMergeEntriesCut; // Should we write the static member?
+  } else {
+    // R__b.WriteClassBuffer(AliTPCcalibAlign::Class(),this);
+    if (isDebug) AliSysInfo::AddStamp("aliengStreamer::Start");
+
+    R__b.WriteVersion(AliTPCcalibAlign::IsA());
+    AliTPCcalibBase::Streamer(R__b); // Stream the base class
+
+    // Write the "big data members directly in the file in parallel with the object itself
+    for (Int_t i=0; i<2; ++i){
+      if (fClusterDelta[i]) fClusterDelta[i]->Write(TString::Format("AliTPCcalibAlign.%s.fClusterDelta_%d",GetName(),i).Data());
+    }
+
+    if (isDebug) AliSysInfo::AddStamp("alignStreamer::fClusterDelta");
+
+    for (Int_t i=0; i<4; ++i){
+      if (fTrackletDelta[i]) fTrackletDelta[i]->Write(TString::Format("AliTPCcalibAlign.%s.fTrackletDelta_%d",GetName(),i).Data());
+    }
+
+    if (isDebug) AliSysInfo::AddStamp("alignStreamer::fTrackletDelta");
+
+    // Write all the other data members. This is error prone,
+    // please make sure this part is updated if you change the data members
+    fDphiHistArray.Streamer(R__b);
+    fDthetaHistArray.Streamer(R__b);
+    fDyHistArray.Streamer(R__b);
+    fDzHistArray.Streamer(R__b);
+
+    fDyPhiHistArray.Streamer(R__b);
+    fDzThetaHistArray.Streamer(R__b);
+    
+    fDphiZHistArray.Streamer(R__b);
+    fDthetaZHistArray.Streamer(R__b);
+    fDyZHistArray.Streamer(R__b);
+    fDzZHistArray.Streamer(R__b);
+      
+    fFitterArray12.Streamer(R__b);
+    fFitterArray9.Streamer(R__b);
+    fFitterArray6.Streamer(R__b);
+
+    fMatrixArray12.Streamer(R__b);
+    fMatrixArray9.Streamer(R__b);
+    fMatrixArray6.Streamer(R__b);
+
+    fCombinedMatrixArray6.Streamer(R__b);
+
+    R__b.WriteFastArray(fPoints,72*72);
+
+    R__b << fNoField;
+    R__b << fXIO;
+    R__b << fXmiddle;
+    R__b << fXquadrant;
+    
+    fArraySectorIntParam.Streamer(R__b);
+    fArraySectorIntCovar.Streamer(R__b);
+
+    R__b << fSectorParamA;
+    R__b << fSectorCovarA;
+    R__b << fSectorParamC;
+    R__b << fSectorCovarC;
+    
+    R__b << fUseInnerOuter;
+    R__b << fgkMergeEntriesCut; // Should we write the static member?
+    if (isDebug) AliSysInfo::AddStamp("alignStreamer::DefaultStreamer");
+  }
+}
index 9cf738cad3cbcf743e2827da5e0556dfca795afb..c48d1bf97bf9791478545119193b81b3ae6908ef 100644 (file)
@@ -206,7 +206,9 @@ protected:
 private:
   AliTPCcalibAlign&  operator=(const AliTPCcalibAlign&);// not implemented
 
-  ClassDef(AliTPCcalibAlign,6)
+  // IMPORTANT: If you change the data members, 
+  // please do not forget to increment the ClassDef and to update the Streamer in AliTPCcalibAlign.cxx
+  ClassDef(AliTPCcalibAlign,7)
 };
 
 
index 584e2ab7878dcc9fe8fd514c6ee4a5b07fb1ae3c..afe28e452b4f0de01389ba37adc8b0773b360752 100644 (file)
@@ -22,7 +22,7 @@
                                                    //   Duplicted also (most probably) AliTPCCalPad --- to be cleaned up
 
 #pragma link C++ class AliTPCcalibCalib+;          // Re-applying calib on cluster level - refitting of the tracks
-#pragma link C++ class AliTPCcalibAlign+;          // (Sector)-alignment calibration 
+#pragma link C++ class AliTPCcalibAlign-;          // (Sector)-alignment calibration 
                                                    //   Histogram cluster to track residuals in order to create later on distortion maps
                                                    // --- update documentation (current documentation is obsolete)