]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Corrections to the custom Streamer
authorfca <fca@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 28 Jun 2000 14:41:12 +0000 (14:41 +0000)
committerfca <fca@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 28 Jun 2000 14:41:12 +0000 (14:41 +0000)
ITS/AliITS.cxx
ITS/AliITS.h

index dce5eebc4038a06c6a696038b5bddbdab920c3ae..724350fe3b908ce56feca5d0f0f9413a5a3acc1c 100644 (file)
@@ -15,6 +15,9 @@
 
 /*
 $Log$
+Revision 1.14  2000/06/15 09:27:52  barbera
+Problems with the HP compiler fixed
+
 Revision 1.13  2000/06/13 15:32:44  nilsen
 fix compilation error on HP and DEC unix.
 
@@ -197,6 +200,12 @@ AliITS::AliITS(const char *name, const char *title):AliDetector(name,title){
    }
   //
 
+   fNdtype = new Int_t[fNDetTypes];
+   fDtype = new TObjArray(fNDetTypes);
+
+   fNctype = new Int_t[fNDetTypes];
+   fCtype = new TObjArray(fNDetTypes);
+
   SetMarkerColor(kRed);
 
   fITSgeom=0;
@@ -668,12 +677,6 @@ void AliITS::MakeBranch(Option_t* option){
   
    char *det[3] = {"SPD","SDD","SSD"};
 
-   fNdtype = new Int_t[fNDetTypes];
-   fDtype = new TObjArray(fNDetTypes);
-
-   fNctype = new Int_t[fNDetTypes];
-   fCtype = new TObjArray(fNDetTypes);
-
    char *kDigclass;
    char *kClclass;
 
@@ -1139,98 +1142,75 @@ Option_t *option,Option_t *opt,Text_t *filename)
 
 }
 
-//____________________________________________________________________________
-void AliITS::Streamer(TBuffer &R__b){
+//______________________________________________________________________________
+void AliITS::Streamer(TBuffer &R__b)
+{
    // Stream an object of class AliITS.
-    Int_t i,j,l;
-
-    AliITSDetType          *det;
-    AliITSresponse         *response;
-    AliITSsegmentation     *seg;
-    TClonesArray           *digitsaddress;
-    TClonesArray           *clustaddr;
 
+  Int_t i, j, l;
 
    if (R__b.IsReading()) {
-      Version_t R__v = R__b.ReadVersion(); 
-      if (R__v == 1) {
-         AliDetector::Streamer(R__b);
-         R__b >> fITSgeom;
-         R__b >> fEuclidOut;
-         R__b >> fIdN;
-         if(fIdSens!=0) delete[] fIdSens;
-         if(fIdName!=0) delete[] fIdName;
-         fIdSens = new Int_t[fIdN];
-         fIdName = new char*[fIdN];
-         for(i=0;i<fIdN;i++) R__b >> fIdSens[i];
-         for(i=0;i<fIdN;i++){
-             R__b >> l;
-             fIdName[i] = new char[l+1]; // add room for null character.
-             for(j=0;j<l;j++) R__b >> fIdName[i][j];
-             fIdName[i][l] = '\0'; // Null terminate this string.
-         } // end for i
-         R__b >> fMajorVersion;
-         R__b >> fMinorVersion;
-         R__b >> fDtype;
-         R__b.ReadArray(fNdtype);
-         R__b >> fCtype;
-         R__b.ReadArray(fNctype);
-         R__b >> fDetTypes;
-         R__b >> fNDetTypes;
-         R__b >> fRecPoints;
-         R__b >> fNRecPoints;
-          //  stream out only response and segmentation
-          for(i=0;i<fNDetTypes;i++) {        
-              det=(AliITSDetType*)(*fDetTypes)[i];
-              det->Streamer(R__b);
-              response=det->GetResponseModel();
-              if (response) response->Streamer(R__b);    
-              seg=det->GetSegmentationModel();
-              if (seg) seg->Streamer(R__b);      
-              digitsaddress=(TClonesArray*) (*fDtype)[i];
-              digitsaddress->Streamer(R__b);
-              clustaddr=(TClonesArray*) (*fCtype)[i];
-              clustaddr->Streamer(R__b);
-         }
-
-      } // end if (R__v)
+      Version_t R__v = R__b.ReadVersion(); if (R__v) { }
+      AliDetector::Streamer(R__b);
+      R__b >> fITSgeom;
+      R__b >> fITSmodules;
+      R__b >> fEuclidOut;
+      R__b >> fIdN;
+      delete []fIdSens; 
+      fIdSens = new Int_t[fIdN]; 
+      R__b.ReadFastArray(fIdSens,fIdN); 
+      for(i=0;i<fIdN;i++){
+       R__b >> l;
+       fIdName[i] = new char[l+1]; // add room for null character.
+       for(j=0;j<l;j++) R__b >> fIdName[i][j];
+       fIdName[i][l] = '\0'; // Null terminate this string.
+      } // end for i
+      //R__b.ReadArray(fIdName);
+      R__b >> fMajorVersion;
+      R__b >> fMinorVersion;
+      R__b >> fDetTypes;
+      R__b >> fNDetTypes;
+      R__b >> fDtype;
+      delete []fNdtype; 
+      fNdtype = new Int_t[fNDetTypes]; 
+      R__b.ReadFastArray(fNdtype,fNDetTypes); 
+      R__b >> fCtype;
+      delete []fNctype; 
+      fNctype = new Int_t[fNDetTypes]; 
+      R__b.ReadFastArray(fNctype,fNDetTypes); 
+      fRecPoints->Streamer(R__b);
+      R__b >> fNRecPoints;
+      fTracks->Streamer(R__b);
+      R__b >> fTreeC;
    } else {
       R__b.WriteVersion(AliITS::IsA());
       AliDetector::Streamer(R__b);
       R__b << fITSgeom;
+      R__b << fITSmodules;
       R__b << fEuclidOut;
       R__b << fIdN;
-      for(i=0;i<fIdN;i++) R__b <<fIdSens[i];
+      R__b.WriteFastArray(fIdSens,fIdN); 
       for(i=0;i<fIdN;i++){
          l = strlen(fIdName[i]);
          R__b << l;
          for(j=0;j<l;j++) R__b << fIdName[i][j];
       } // end for i
+      //R__b.WriteArray(fIdName, __COUNTER__);
       R__b << fMajorVersion;
       R__b << fMinorVersion;
-      R__b << fDtype;
-      R__b.WriteArray(fNdtype, fNDetTypes);
-      R__b << fCtype;
-      R__b.WriteArray(fNctype, fNDetTypes);
       R__b << fDetTypes;
       R__b << fNDetTypes;
+      R__b << fDtype;
+      R__b.WriteFastArray(fNdtype,fNDetTypes); 
+      R__b << fCtype;
+      R__b.WriteFastArray(fNctype,fNDetTypes); 
       R__b << fRecPoints;
+      // fRecPoints->Streamer(R__b);
       R__b << fNRecPoints;
-      for(i=0;i<fNDetTypes;i++) {
-          det=(AliITSDetType*)(*fDetTypes)[i];
-          det->Streamer(R__b);
-           response=det->GetResponseModel();
-          if(response) response->Streamer(R__b);
-          seg=det->GetSegmentationModel();
-          if (seg) seg->Streamer(R__b);          
-          digitsaddress=(TClonesArray*) (*fDtype)[i];
-          digitsaddress->Streamer(R__b);
-          clustaddr=(TClonesArray*) (*fCtype)[i];
-          clustaddr->Streamer(R__b);
-      }          
+      R__b << fTracks;
+      //      fTracks->Streamer(R__b);
+      R__b << fTreeC;
    }
-
-
 }
 
 ClassImp(AliITSRecPoint)
index c581b4f5c5c6a46bff03a9147e40d19246c320d8..cfc4bd50b0fcbf1a065159b89a8937e18583adfb 100644 (file)
@@ -145,7 +145,8 @@ class AliITS : public AliDetector {
     Bool_t fEuclidOut;                       // Flag to write out geometry 
                                              //in euclid format
     Int_t  fIdN;                             // the number of layers
-    Int_t  *fIdSens; char **fIdName;         //layer identifier
+    Int_t  *fIdSens;                         //[fIdN] layer identifier
+    char **fIdName;                          //layer identifier
     // Geometry and Stepmanager version numbers used.
     Int_t fMajorVersion,fMinorVersion;       //detailed and coarse(minor) versions
     //
@@ -154,10 +155,10 @@ class AliITS : public AliDetector {
     Int_t                 fNDetTypes;          // Number of Detector types
 
     TObjArray            *fDtype;              // List of digits
-    Int_t                *fNdtype;             // Number of digits per type of
+    Int_t                *fNdtype;             //[fNDetTypes] Number of digits per type of
                                                // detector 
     TObjArray            *fCtype;              // List of clusters
-    Int_t                *fNctype;             // Number of clusters per type
+    Int_t                *fNctype;             //[fNDetTypes] Number of clusters per type
                                                // of detector