]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - ITS/AliITSDDLRawData.cxx
load AddTaskEmcalJet.C
[u/mrichter/AliRoot.git] / ITS / AliITSDDLRawData.cxx
index 79e11cc8b9b305b77b8f2a7228038152f3d8df47..a0ba53ead55bc81492bb327d6cc75945f643d178 100644 (file)
@@ -39,6 +39,9 @@
 #include "AliFstream.h"
 #include "AliITSFOSignalsSPD.h"
 
+using std::ofstream;
+using std::ios;
+using std::endl;
 ClassImp(AliITSDDLRawData)
 
 ////////////////////////////////////////////////////////////////////////////////////////
@@ -46,7 +49,7 @@ AliITSDDLRawData::AliITSDDLRawData():
 fVerbose(0),
 fIndex(-1),
 fHalfStaveModule(-1),
-fUseCompressedSDDFormat(0){
+fSDDRawFormat(7){
   //Default constructor
 
 }
@@ -58,7 +61,7 @@ AliITSDDLRawData::AliITSDDLRawData(const AliITSDDLRawData &source) :
 fVerbose(source.fVerbose),
 fIndex(source.fIndex),
 fHalfStaveModule(source.fHalfStaveModule),
-fUseCompressedSDDFormat(source.fUseCompressedSDDFormat){
+fSDDRawFormat(source.fSDDRawFormat){
   //Copy Constructor
 }
 
@@ -66,10 +69,11 @@ fUseCompressedSDDFormat(source.fUseCompressedSDDFormat){
 
 AliITSDDLRawData& AliITSDDLRawData::operator=(const AliITSDDLRawData &source){
   //Assigment operator
-  this->fIndex=source.fIndex;
-  this->fHalfStaveModule=source.fHalfStaveModule;
-  this->fVerbose=source.fVerbose;
-  this->fUseCompressedSDDFormat=source.fUseCompressedSDDFormat;
+  if(this==&source) return *this;
+  fIndex=source.fIndex;
+  fHalfStaveModule=source.fHalfStaveModule;
+  fVerbose=source.fVerbose;
+  fSDDRawFormat=source.fSDDRawFormat;
   return *this;
 }
 
@@ -98,7 +102,8 @@ void AliITSDDLRawData::GetDigitsSSD(TClonesArray *ITSdigits,Int_t mod,Int_t modR
       ix=digs->GetCoord2();  // Strip Number
       is=digs->GetCompressedSignal();  // ADC Signal
       // cout<<" Module:"<<mod-500<<" N/P side:"<<iz<<" Strip Number:"<<ix<<" Amplidute:"<<is-1<<endl;
-      if(is<0) is = 4096 + is;
+      if(is<0) is = 0;
+      if(is>4095) is = 4095;
       if (fVerbose==2)
        ftxt<<"DDL:"<<ddl<<" Mod: "<<modR<<" N/P: "<<iz<<" Strip: "<<ix<<" Value: "<<is-1<<endl;
 
@@ -151,12 +156,21 @@ void AliITSDDLRawData::GetDigitsSDDCompressed(TClonesArray *ITSdigits, Int_t mod
       dataWord+=sid<<26;
       dataWord+=iz<<18;
       dataWord+=ix<<10;
-      dataWord+=is;
+      UInt_t adcEncoded=0;
+      Int_t shift=0;
+      if(is < 8) shift=2;
+      else if(is<16) shift=3;
+      else if(is<32) shift=4;
+      else if(is<64) shift=5;
+      else if(is<128) shift=6;
+      else shift=7;
+      adcEncoded=shift+((is-(1<<shift))<<3);
+      dataWord+=adcEncoded;
       fIndex++;
       buf[fIndex]=dataWord;
     }
   }
-  UInt_t finalWord=15<<28;
+  UInt_t finalWord=UInt_t(15)<<28;
   finalWord+=mod;
   fIndex++;
   buf[fIndex]=finalWord;  
@@ -534,14 +548,14 @@ Int_t AliITSDDLRawData::RawDataSPD(TBranch* branch, AliITSFOSignalsSPD* foSignal
   fIndex=-1;
 
   TClonesArray*& digits = * (TClonesArray**) branch->GetAddress();
-  char fileName[15];
+  TString fileName;
   AliFstream* outfile;         // logical name of the output file 
   AliRawDataHeaderSim header;
 
   //loop over DDLs
   for(Int_t ddl=0;ddl<AliDAQ::NumberOfDdls("ITSSPD");ddl++){
-    strcpy(fileName,AliDAQ::DdlFileName("ITSSPD",ddl)); //The name of the output file.
-    outfile = new AliFstream(fileName);
+    fileName.Form("%s",AliDAQ::DdlFileName("ITSSPD",ddl)); //The name of the output file.
+    outfile = new AliFstream(fileName.Data());
     //write Dummy DATA HEADER
     UInt_t dataHeaderPosition=outfile->Tellp();
     outfile->WriteBuffer((char*)(&header),sizeof(header));
@@ -581,14 +595,14 @@ Int_t AliITSDDLRawData::RawDataSSD(TBranch* branch){
   fIndex=-1;
 
   TClonesArray*& digits = * (TClonesArray**) branch->GetAddress();
-  char fileName[15];
+  TString fileName;
   AliFstream* outfile;         // logical name of the output file 
   AliRawDataHeaderSim header;
 
   //loop over DDLs  
   for(Int_t i=0;i<AliDAQ::NumberOfDdls("ITSSSD");i++){
-    strcpy(fileName,AliDAQ::DdlFileName("ITSSSD",i)); //The name of the output file.
-    outfile = new AliFstream(fileName);
+    fileName.Form("%s",AliDAQ::DdlFileName("ITSSSD",i)); //The name of the output file.
+    outfile = new AliFstream(fileName.Data());
     //write Dummy DATA HEADER
     UInt_t dataHeaderPosition=outfile->Tellp();
     outfile->WriteBuffer((char*)(&header),sizeof(header));
@@ -621,34 +635,45 @@ Int_t AliITSDDLRawData::RawDataSSD(TBranch* branch){
 
 /////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
 
-Int_t AliITSDDLRawData::RawDataSDD(TBranch* branch, AliITSDDLModuleMapSDD* ddlsdd){
+Int_t AliITSDDLRawData::RawDataSDD(TBranch* branch, const AliITSDDLModuleMapSDD* ddlsdd){
     //This method creates the Raw data files for SDD detectors
   const Int_t kSize=131072; //256*512
   UInt_t buf[kSize];      
   fIndex=-1;
 
   TClonesArray*& digits = * (TClonesArray**) branch->GetAddress();
-  char fileName[15];
+  TString fileName;
   AliFstream* outfile;             // logical name of the output file 
   AliRawDataHeaderSim header;
-  UInt_t skippedword, carlosFooterWord,fifoFooterWord,jitterWord;
+  
+  if(fSDDRawFormat!=0){ 
+    for(Int_t ibit=0; ibit<8; ibit++) header.SetAttribute(ibit);
+  }else{
+    for(Int_t ibit=0; ibit<5; ibit++) header.SetAttribute(ibit);
+    for(Int_t ibit=5; ibit<8; ibit++) header.ResetAttribute(ibit);  
+  }
+  UInt_t skippedword=0; 
+  UInt_t carlosFooterWord=0;
+  UInt_t fifoFooterWord=0;
+  UInt_t jitterWord=0;
   Bool_t retcode;
   retcode = AliBitPacking::PackWord(0x3FFFFFFF,carlosFooterWord,0,31);
+  if(!retcode)AliError("AliBitPacking error\n");
   retcode = AliBitPacking::PackWord(0x3F1F1F1F,fifoFooterWord,0,31);
-  if(!fUseCompressedSDDFormat) retcode = AliBitPacking::PackWord(0x7F000000,jitterWord,0,31);
+  if(fSDDRawFormat!=0) retcode = AliBitPacking::PackWord(0x7F000000,jitterWord,0,31);
   else retcode = AliBitPacking::PackWord(0x80000000,jitterWord,0,31);
  
   //loop over DDLs  
   for(Int_t i=0;i<AliDAQ::NumberOfDdls("ITSSDD");i++){
-    strcpy(fileName,AliDAQ::DdlFileName("ITSSDD",i)); //The name of the output file.
-    outfile = new AliFstream(fileName);
+    fileName.Form("%s",AliDAQ::DdlFileName("ITSSDD",i)); //The name of the output file.
+    outfile = new AliFstream(fileName.Data());
     //write Dummy DATA HEADER
     UInt_t dataHeaderPosition=outfile->Tellp();
     outfile->WriteBuffer((char*)(&header),sizeof(header));
 
 
     //first 1 "dummy" word to be skipped
-    if(!fUseCompressedSDDFormat){
+    if(fSDDRawFormat!=0){
       retcode = AliBitPacking::PackWord(0xFFFFFFFF,skippedword,0,31);
       outfile->WriteBuffer((char*)(&skippedword),sizeof(skippedword));
     }
@@ -663,7 +688,7 @@ Int_t AliITSDDLRawData::RawDataSDD(TBranch* branch, AliITSDDLModuleMapSDD* ddlsd
        //For each Module, buf contains the array of data words in Binary format          
        //fIndex gives the number of 32 bits words in the buffer for each module
        //      cout<<"MODULE NUMBER:"<<mapSDD[i][mod]<<endl;
-       if(fUseCompressedSDDFormat){
+       if(fSDDRawFormat==0){
          GetDigitsSDDCompressed(digits,mod,buf);
          outfile->WriteBuffer((char *)buf,((fIndex+1)*sizeof(UInt_t)));
        }else{
@@ -675,7 +700,7 @@ Int_t AliITSDDLRawData::RawDataSDD(TBranch* branch, AliITSDDLModuleMapSDD* ddlsd
       }//end if
     }//end for
     // 12 words with FIFO footers (=4 FIFO x 3 3F1F1F1F words per DDL)
-    if(!fUseCompressedSDDFormat){
+    if(fSDDRawFormat!=0){
       for(Int_t iw=0;iw<12;iw++) outfile->WriteBuffer((char*)(&fifoFooterWord),sizeof(fifoFooterWord));
     }
     outfile->WriteBuffer((char*)(&jitterWord),sizeof(jitterWord));      
@@ -742,3 +767,4 @@ void  AliITSDDLRawData::WriteHit(UInt_t *buf,Int_t RowAddr,Int_t HitAddr,UInt_t
   }//end else
   return;
 }//end WriteHit
+