]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
SDD Quality Assurance protoype macro - SDD RawStreamer code changed in order to be...
authormasera <masera@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sat, 25 Aug 2007 16:18:25 +0000 (16:18 +0000)
committermasera <masera@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sat, 25 Aug 2007 16:18:25 +0000 (16:18 +0000)
ITS/AliITSRawStreamSDD.cxx
ITS/AliITSRawStreamSDD.h
ITS/ITSSDDQA.C [new file with mode: 0644]

index bdecf5cf1445bb9fdb11f2a9d30e4e7285afaffc..15e3ae56ada363d7860a0f66c1f6d45faeb5fcdb 100644 (file)
@@ -24,6 +24,7 @@
 #include "AliITSRawStreamSDD.h"
 #include "AliRawReader.h"
 #include "AliLog.h"
+#include "Riostream.h"
 
 ClassImp(AliITSRawStreamSDD)
 
@@ -64,10 +65,6 @@ fEventId(0),
 fChannel(0),
 fJitter(0),
 fNCarlos(kModulesPerDDL),
-fNfifo0(0),
-fNfifo1(0),
-fNfifo2(0),
-fNfifo3(0),
 fDDL(0){
 // create an object to read ITS SDD raw digits
   
@@ -83,12 +80,14 @@ fDDL(0){
     }
     fLowThreshold[i]=0;
   }
+  for(Int_t i=0;i<kFifoWords;i++) fNfifo[i]=0;
   for(Int_t i=0;i<kDDLsNumber;i++) fSkip[i]=0;
   fRawReader->Reset();
   fRawReader->Select("ITSSDD");
 
-  //  fRawReader->SelectEquipment(17, 101, 101);//select this for test data
   //  fNCarlos = 8; //select this for test data
+  for(Short_t i=0; i<kCarlosWords; i++) iCarlosWord[i]=805306368+i;
+  for(Short_t i=0; i<kFifoWords; i++) iFifoWord[i]=805306384+i;
 }
 
 UInt_t AliITSRawStreamSDD::ReadBits()
@@ -125,31 +124,11 @@ Bool_t AliITSRawStreamSDD::Next()
   fDDL=fRawReader->GetDDLID();
   Int_t ddln = fRawReader->GetDDLID();
   if(ddln <0) ddln=0;
-  while (fSkip[ddln] < 9) {
-    if (!fRawReader->ReadNextInt(fData))return kFALSE;
-    if ((fData >> 30) == 0x01) continue;  // JTAG word
-    fSkip[ddln]++;
-  }
-  
+  Bool_t kSkip = ResetSkip(ddln);
+  if(!kSkip) return kSkip;
+
+  for(Int_t i=0;i<kModulesPerDDL;i++) iCountFoot[i]=0; 
 
-  Int_t countFoot[kModulesPerDDL];
-  for(Int_t i=0;i<kModulesPerDDL;i++) countFoot[i]=0;
-  UInt_t iCarlos0Word=805306368;
-  UInt_t iCarlos1Word=805306369;
-  UInt_t iCarlos2Word=805306370;
-  UInt_t iCarlos3Word=805306371;
-  UInt_t iCarlos4Word=805306372;
-  UInt_t iCarlos5Word=805306373;
-  UInt_t iCarlos6Word=805306374;
-  UInt_t iCarlos7Word=805306375;
-  UInt_t iCarlos8Word=805306376;
-  UInt_t iCarlos9Word=805306377;
-  UInt_t iCarlos10Word=805306378;
-  UInt_t iCarlos11Word=805306379;
-  UInt_t iFifo0Word=805306384;
-  UInt_t iFifo1Word=805306385;
-  UInt_t iFifo2Word=805306386;
-  UInt_t iFifo3Word=805306387;
   while (kTRUE) {
     if ((fChannel < 0) || (fLastBit[fCarlosId][fChannel] < fReadBits[fCarlosId][fChannel])) {
 
@@ -157,86 +136,87 @@ Bool_t AliITSRawStreamSDD::Next()
        return kFALSE;  // read next word 
       }
       ddln = fRawReader->GetDDLID();
-      if(ddln!=fDDL){
-       Reset();
-       for(Int_t icr=0;icr<kModulesPerDDL;icr++) countFoot[icr]=0;
-      }
+
+      if(ddln!=fDDL) Reset();
       if(ddln < 0 || ddln > (kDDLsNumber-1)) ddln  = 0;
-      while (fSkip[ddln] < 9) {
-       if (!fRawReader->ReadNextInt(fData)){ 
-         return kFALSE;
-       }
-       if ((fData >> 30) == 0x01) continue;  // JTAG word
-       fSkip[ddln]++;
-      }
+      kSkip = ResetSkip(ddln);
+      if(!kSkip) return kFALSE;  // read next word
 
       fChannel = -1;
-      if((fData>=iCarlos0Word&&fData<=iCarlos11Word)||(fData>=iFifo0Word&&fData<=iFifo3Word)){
-       if(fData==iCarlos0Word){
+      if((fData>=iCarlosWord[0]&&fData<=iCarlosWord[11])||(fData>=iFifoWord[0]&&fData<=iFifoWord[3])){
+       for(Short_t ik=0;ik<kCarlosWords;ik++) {
+         if(fData==iCarlosWord[ik]) {
+           fCarlosId = ik;
+           Int_t iFifoIdx = ik/3;
+           if(fNCarlos == 8) {
+             if(ik==2) iFifoIdx = 1;
+             if(ik==4 || ik==5)  iFifoIdx = 2;
+             if(ik==6 || ik==7)  iFifoIdx = 3 ;
+           }       
+           fNfifo[iFifoIdx] = fCarlosId;
+         }
+       }
+       for(Short_t ik=0;ik<kFifoWords;ik++) {
+         if (fData==iFifoWord[ik]) { 
+           fCarlosId = fNfifo[ik];         
+         }
+       }
+       /*
+       if(fData==iCarlosWord[0]){
          fCarlosId = 0;
-         fNfifo0 = fCarlosId;
+         fNfifo[0] = fCarlosId;
        }
-       else if (fData==iCarlos1Word){
+       else if (fData==iCarlosWord[1]){
          fCarlosId = 1;
-         fNfifo0 = fCarlosId;
+         fNfifo[0] = fCarlosId;
        }
-       else if (fData==iCarlos2Word){
+       else if (fData==iCarlosWord[2]){
          fCarlosId = 2;
-         if(fNCarlos == 8) fNfifo1 = fCarlosId;
-         else fNfifo0 = fCarlosId;
+         if(fNCarlos == 8) fNfifo[1] = fCarlosId;
+         else fNfifo[0] = fCarlosId;
        }
-       else if (fData==iCarlos3Word){
+       else if (fData==iCarlosWord[3]){
          fCarlosId = 3;
-         fNfifo1 = fCarlosId;
+         fNfifo[1] = fCarlosId;
        }
-       else if (fData==iCarlos4Word){
+       else if (fData==iCarlosWord[4]){
          fCarlosId = 4;
-         if(fNCarlos == 8) fNfifo2 = fCarlosId;
-         else fNfifo1 = fCarlosId;
+         if(fNCarlos == 8) fNfifo[2] = fCarlosId;
+         else fNfifo[1] = fCarlosId;
        }
-       else if (fData==iCarlos5Word){
+       else if (fData==iCarlosWord[5]){
          fCarlosId = 5;
-         if(fNCarlos == 8) fNfifo2 = fCarlosId;
-         else fNfifo1 = fCarlosId;
+         if(fNCarlos == 8) fNfifo[2] = fCarlosId;
+         else fNfifo[1] = fCarlosId;
        }
-       else if (fData==iCarlos6Word){
+       else if (fData==iCarlosWord[6]){
          fCarlosId = 6;
-         if(fNCarlos == 8) fNfifo3 = fCarlosId;
-         else fNfifo2 = fCarlosId;
+         if(fNCarlos == 8) fNfifo[3] = fCarlosId;
+         else fNfifo[2] = fCarlosId;
        }
-       else if (fData==iCarlos7Word){
+       else if (fData==iCarlosWord[7]){
          fCarlosId = 7;
-         if(fNCarlos == 8) fNfifo3 = fCarlosId;
-         else fNfifo2 = fCarlosId;
+         if(fNCarlos == 8) fNfifo[3] = fCarlosId;
+         else fNfifo[2] = fCarlosId;
        }
-       else if (fData==iCarlos8Word){
+       else if (fData==iCarlosWord[8]){
          fCarlosId = 8;
-         fNfifo2 = fCarlosId;
+         fNfifo[2] = fCarlosId;
        }
-       else if (fData==iCarlos9Word){
+       else if (fData==iCarlosWord[9]){
          fCarlosId = 9;
-         fNfifo3 = fCarlosId;
+         fNfifo[3] = fCarlosId;
        }
-       else if (fData==iCarlos10Word){
+       else if (fData==iCarlosWord[10]){
          fCarlosId = 10;
-         fNfifo3 = fCarlosId;
+         fNfifo[3] = fCarlosId;
        }
-       else if (fData==iCarlos11Word){
+       else if (fData==iCarlosWord[11]){
          fCarlosId = 11;
-         fNfifo3 = fCarlosId;
-       }
-       else if (fData==iFifo0Word){
-         fCarlosId = fNfifo0;
-       }
-       else if (fData==iFifo1Word){
-         fCarlosId = fNfifo1;
-       }
-       else if (fData==iFifo2Word){
-         fCarlosId = fNfifo2;
-       }
-       else if (fData==iFifo3Word){
-         fCarlosId = fNfifo3;
+         fNfifo[3] = fCarlosId;
        }
+       */
+
       }
       if(fData==1059004191) continue;
       if (fNCarlos == 8 && (fCarlosId == 8 || fCarlosId == 9 || 
@@ -247,10 +227,10 @@ Bool_t AliITSRawStreamSDD::Next()
       if ((fData >> 28) == 0x02) {           // header
        fEventId = (fData >> 3) & 0x07FF;
       } else if ((fData >> 28) == 0x03) {    // footer
-        countFoot[fCarlosId]++; // stop before the last word (last word=jitter)
+        iCountFoot[fCarlosId]++; // stop before the last word (last word=jitter)
        Bool_t exit=kTRUE;
        for(Int_t ic=0;ic<fNCarlos;ic++){
-         if(countFoot[ic]<3) exit=kFALSE;
+         if(iCountFoot[ic]<3) exit=kFALSE;
        }
         if(exit){
          return kFALSE;
@@ -313,8 +293,8 @@ Bool_t AliITSRawStreamSDD::Next()
 }
 
 void AliITSRawStreamSDD::Reset(){
+
   //reset data member for a new ddl
-  
   for(Int_t i=0;i<2;i++){
     for(Int_t ic=0;ic<kModulesPerDDL;ic++){
       fChannelData[ic][i]=0;
@@ -327,6 +307,18 @@ void AliITSRawStreamSDD::Reset(){
     }
     fLowThreshold[i]=0;
   }
+  for(Int_t i=0;i<kModulesPerDDL;i++) iCountFoot[i]=0;
+  
   fChannel=-1;
   fDDL=fRawReader->GetDDLID();
 }
+
+Bool_t AliITSRawStreamSDD::ResetSkip(Int_t ddln){
+  while (fSkip[ddln] < 9) {
+    if (!fRawReader->ReadNextInt(fData))return kFALSE;
+    if ((fData >> 30) == 0x01) continue;  // JTAG word
+    fSkip[ddln]++;
+  }
+  return kTRUE;
+}
+
index b57989c21b57db0a4280e2bfdc5dbf84eb12aeeb..72a0fce6fcb69a80055071d05f6a2703e04bfb43 100644 (file)
@@ -33,10 +33,13 @@ class AliITSRawStreamSDD: public AliITSRawStream {
     static  Int_t    GetModuleNumber(UInt_t iDDL, UInt_t iModule)
                      {return fgkDDLModuleMap[iDDL][iModule];}
     virtual void     Reset(); 
+    virtual Bool_t   ResetSkip(Int_t ddln); 
 
     enum {kDDLsNumber = 24};      // number of DDLs in SDD
     enum {kModulesPerDDL = 12};   // number of modules in each DDL 
-    enum ESDDRawStreamError {
+    enum {kCarlosWords = 12};      // number of FIFOCARLOS Words
+    enum {kFifoWords =  4};      // number of FIFO Words
+    enum ESDDRawStreamError { 
       kDataError = 1,
       kDataFormatErr = 2
     };
@@ -61,15 +64,20 @@ class AliITSRawStreamSDD: public AliITSRawStream {
     UInt_t           fReadBits[kModulesPerDDL][2];   // number of bits to read
     Int_t            fLowThreshold[2]; // low Carlos threshold
     Int_t            fNCarlos;         // number of Carlos 
+    /*
     Int_t            fNfifo0;          // fifo n. 0
     Int_t            fNfifo1;          // fifo n. 1
     Int_t            fNfifo2;          // fifo n. 2
     Int_t            fNfifo3;          // fifo n. 3
+    */
+    Int_t            fNfifo[kFifoWords];
     Int_t            fTimeBin[kModulesPerDDL][2];  // current timebin [ncarlos][nchannels]
     Int_t            fAnode[kModulesPerDDL][2]; // current anode [ncarlos][nchannels]
     Int_t            fDDL;        //current ddl number
-
-    ClassDef(AliITSRawStreamSDD, 2) // class for reading ITS SDD raw digits
+    UInt_t           iCarlosWord[kCarlosWords];
+    UInt_t           iFifoWord[kFifoWords];
+    Int_t            iCountFoot[kModulesPerDDL];
+    ClassDef(AliITSRawStreamSDD, 4) // class for reading ITS SDD raw digits
 };
 
 #endif
diff --git a/ITS/ITSSDDQA.C b/ITS/ITSSDDQA.C
new file mode 100644 (file)
index 0000000..cd8ed45
--- /dev/null
@@ -0,0 +1,126 @@
+void ITSSDDQA(char *iFile, Int_t MaxEvts=1000000, Int_t FirstEvt=0) {
+
+cout << "SDD Quality Assurance Prototype Macro" << endl; 
+
+const Int_t nSDDmodules= 260;
+const Int_t imodoffset = 240;
+const Int_t modtotSDD  = nSDDmodules*2;
+
+Float_t xi = 0.5;
+Float_t xf = xi + nSDDmodules;
+TH1F *modulePattern = new TH1F("patternModule","Modules pattern",nSDDmodules,xi,xf);
+xf = xi + modtotSDD;
+TH1F *moduleSidePattern = new TH1F("patternSide","Modules/Side pattern",modtotSDD,xi,xf);
+
+TH2F *hismap[modtotSDD];  //260 dx e 260 sx  with A, T, Q
+TH2F *hispop[modtotSDD];  //260 dx e 260 sx  with A, T, Ncounts
+Char_t *cindex = new Char_t[5];
+for(Int_t imod=0; imod<nSDDmodules;imod++){
+  for(Int_t isid=0;isid<2;isid++){
+    Int_t index=2*imod+isid;       //260*2 position
+
+    sprintf(cindex,"%d",index+1); // imod,isid);
+    TString sindex((const char *) cindex);
+
+    TString histnam("chargeMap");
+    TString histit("Total Charge, module number ");
+    histnam.Append(sindex);
+    histit.Append(sindex);
+    hismap[index]=new TH2F(histnam.Data(),histit.Data(),256,-0.5,255.5,256,-0.5,255.5);
+      
+    TString hisnam2("countsMap");
+    TString histit2("Number of Counts, module number ");
+    hisnam2.Append(sindex);
+    histit2.Append(sindex);
+    hispop[index]=new TH2F(hisnam2.Data(),histit2.Data(),256,-0.5,255.5,256,-0.5,255.5);
+  
+    /*
+      sprintf(hisnam,"hisprojX%03ds%d",imod,isid);
+      sprintf(histitle,"layer , ladder, module position, channel, %d, %d", imod, isid);
+      hisprojX[index]=new TH2F(hisnam,histitle,256,-0.5,255.5,256,-0.5,255.5);
+
+      sprintf(hisnam,"hisprojY%03ds%d",imod,isid);
+      sprintf(histitle,"layer , ladder, module position, channel, %d, %d", imod, isid);
+      hisprojY[index]=new TH2F(hisnam,histitle,256,-0.5,255.5,256,-0.5,255.5);
+
+      sprintf(hisnam,"hisprofX%03ds%d",imod,isid);
+      sprintf(histitle,"layer , ladder, module position, channel, %d, %d", imod, isid);
+      hisprofX[index]=new TH2F(hisnam,histitle,256,-0.5,255.5,256,-0.5,255.5);
+
+      sprintf(hisnam,"hisprofY%03ds%d",imod,isid);
+      sprintf(histitle,"layer , ladder, module position, channel, %d, %d", imod, isid);
+      hisprofY[index]=new TH2F(hisnam,histitle,256,-0.5,255.5,256,-0.5,255.5);
+      */
+  }
+}
+delete [] cindex;
+
+  AliRawReader *rd = new AliRawReaderDate(iFile,FirstEvt);  // open run
+  Int_t evCounter = 0;
+
+  //AliITS *itsRun = new AliITS();
+  //  TGeoManager::Import("$ALICE_ROOT/EVE/alice-data/alice_fullgeo.root");
+  /*
+    AliITSInitGeometry *initgeom = new AliITSInitGeometry("AliITSvPPRasymmFMD",2);
+    geom = initgeom->CreateAliITSgeom();
+    delete initgeom;
+  */
+
+  Int_t eqOffset = 256;
+  Int_t DDLid_range = 24;
+  do{       // start loop on events
+    if(++evCounter > MaxEvts) { cout << MaxEvts << " events read, stop" << endl; evCounter--; break; }  
+     cout << "Read Event: " << evCounter+FirstEvt-1 << endl;
+
+    rd->RequireHeader(kFALSE);             
+    rd->SelectEvents(7);                   // read only events with the given type. no selection is applied if a value < 0 is used. 
+
+    rd->SelectEquipment(17,eqOffset+1,eqOffset+DDLid_range);
+
+    rd->Reset();    // reset the current position to the beginning of the event
+    AliITSRawStreamSDD s(rd);    //This class provides access to ITS SDD digits in raw data.
+    Int_t iddl;
+    Int_t isddmod;
+    Int_t moduleSDD;
+    while(s.Next()){                       //read the next raw digit; returns kFALSE if there is no digit left
+
+       iddl=rd->GetDDLID();
+       isddmod=s.GetModuleNumber(iddl,s.GetCarlosId());  //this is the FEE Carlos
+       //cout<<"DDLID= "<<iddl <<"; Module number= " <<isddmod  <<endl;
+       modulePattern->Fill(isddmod-imodoffset+1);     // 1 to 260
+       moduleSDD=2*(isddmod-imodoffset)+s.GetChannel();
+        moduleSidePattern->Fill(moduleSDD+1);          // 1 to 520
+        //cout << "anode " << s.GetCoord1() << ", time bin: " << s.GetCoord2() << ", charge: " << s.GetSignal() << endl;
+        hismap[moduleSDD]->Fill(s.GetCoord1(), s.GetCoord2(),s.GetSignal() );
+        hispop[moduleSDD]->Fill(s.GetCoord1(), s.GetCoord2() );         
+    }    
+
+  } while(rd->NextEvent()); // end loop on events
+  delete rd;
+
+  cout << "end after " << evCounter << " events" << endl;
+
+   TString oFileName(iFile);
+   oFileName.Append(".root");
+   TFile *oFile = TFile::Open(oFileName,"recreate");
+   modulePattern->Write();
+   moduleSidePattern->Write();
+   for(Int_t i=0; i<modtotSDD; i++){
+     hismap[i]->Write();
+     hispop[i]->Write();
+   }
+   oFile->Close();
+
+   for(Int_t imod=0; imod<nSDDmodules;imod++){
+     for(Int_t isid=0;isid<2;isid++){
+       Int_t index=2*imod+isid;       //260*2 position
+       delete hismap[index];
+       delete hispop[index];
+     }
+   }
+   delete modulePattern;
+   delete moduleSidePattern;
+
+}