]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Allow up-to oversampling 4
authorcholm <cholm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 26 Nov 2007 13:50:44 +0000 (13:50 +0000)
committercholm <cholm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 26 Nov 2007 13:50:44 +0000 (13:50 +0000)
21 files changed:
FMD/AliFMD.cxx
FMD/AliFMD.h
FMD/AliFMDBaseDigitizer.cxx
FMD/AliFMDBaseDigitizer.h
FMD/AliFMDCalibFaker.cxx
FMD/AliFMDDigit.cxx
FMD/AliFMDDigit.h
FMD/AliFMDDigitizer.cxx
FMD/AliFMDDigitizer.h
FMD/AliFMDRawReader.cxx
FMD/AliFMDRawStream.cxx
FMD/AliFMDReconstructor.cxx
FMD/AliFMDSDigit.cxx
FMD/AliFMDSDigit.h
FMD/AliFMDSDigitizer.cxx
FMD/AliFMDSDigitizer.h
FMD/scripts/DrawHits.C
FMD/scripts/MakeCalibration.C
FMD/scripts/RawTest.C
FMD/scripts/ReadRaw.C
FMD/scripts/TestShaping.C [new file with mode: 0644]

index 98e5fa9b1b837ee1b90e450c3737f58d35ee5f02..ee59f3846cb00631129902d754bf3c422a9a660e 100644 (file)
@@ -879,7 +879,8 @@ AliFMD::AddDigit(Int_t* digits, Int_t*)
                   UShort_t(digits[3]),  // Strip #
                   UShort_t(digits[4]),  // ADC Count1 
                   Short_t(digits[5]),   // ADC Count2 
-                  Short_t(digits[6]));  // ADC Count3 
+                  Short_t(digits[6]),   // ADC Count3 
+                  Short_t(digits[7])); 
 }
 
 //____________________________________________________________________
@@ -890,7 +891,8 @@ AliFMD::AddDigitByFields(UShort_t detector,
                         UShort_t strip, 
                         UShort_t count1, 
                         Short_t  count2,
-                        Short_t  count3)
+                        Short_t  count3, 
+                        Short_t  count4)
 {
   // add a real digit - as coming from data
   // 
@@ -906,10 +908,11 @@ AliFMD::AddDigitByFields(UShort_t detector,
   TClonesArray& a = *(DigitsArray());
   
   new (a[fNdigits++]) 
-    AliFMDDigit(detector, ring, sector, strip, count1, count2, count3);
-  AliFMDDebug(15, ("Adding digit # %5d/%5d for FMD%d%c[%2d,%3d]=(%d,%d,%d)",
+    AliFMDDigit(detector, ring, sector, strip, count1, count2, count3, count4);
+  AliFMDDebug(15, ("Adding digit # %5d/%5d for FMD%d%c[%2d,%3d]=(%d,%d,%d,%d)",
                   fNdigits-1, a.GetEntriesFast(),
-                  detector, ring, sector, strip, count1, count2, count3));
+                  detector, ring, sector, strip, 
+                  count1, count2, count3, count4));
   
 }
 
@@ -937,7 +940,8 @@ AliFMD::AddSDigit(Int_t* digits)
                    Float_t(digits[4]),   // Edep
                    UShort_t(digits[5]),  // ADC Count1 
                    Short_t(digits[6]),   // ADC Count2 
-                   Short_t(digits[7]));  // ADC Count3 
+                   Short_t(digits[7]),   // ADC Count3 
+                   Short_t(digits[8]));
 }
 
 //____________________________________________________________________
@@ -949,7 +953,8 @@ AliFMD::AddSDigitByFields(UShort_t detector,
                          Float_t  edep,
                          UShort_t count1, 
                          Short_t  count2,
-                         Short_t  count3)
+                         Short_t  count3, 
+                         Short_t  count4)
 {
   // add a summable digit
   // 
@@ -967,7 +972,8 @@ AliFMD::AddSDigitByFields(UShort_t detector,
   TClonesArray& a = *(SDigitsArray());
   
   new (a[fNsdigits++]) 
-    AliFMDSDigit(detector, ring, sector, strip, edep, count1, count2, count3);
+    AliFMDSDigit(detector, ring, sector, strip, edep, 
+                count1, count2, count3, count4);
 }
 
 //____________________________________________________________________
index db18d7e8e013741a0bd51b536fcffdb45de34651..13a611e5badeb8458e28a16f1ff10c16f0eb5c94 100644 (file)
@@ -445,7 +445,8 @@ public:
                                         UShort_t strip=0, 
                                         UShort_t count1=0, 
                                         Short_t  count2=-1, 
-                                        Short_t  count3=-1);
+                                        Short_t  count3=-1, 
+                                        Short_t  count4=-1);
   /** Add a digit to the Digit tree 
       @param digits
       - digits[0]  [UShort_t] Detector #
@@ -472,7 +473,8 @@ public:
                                          Float_t  edep=0,
                                          UShort_t count1=0, 
                                          Short_t  count2=-1, 
-                                         Short_t  count3=-1);
+                                         Short_t  count3=-1,
+                                         Short_t  count4=-1);
   /** @}*/
 
   /** @{ */
index 11c7056b100e458f863c5c320b436af6580b1330..fcb93a3da838465297f5b0b59cebc3e0875ecbe1 100644 (file)
@@ -227,7 +227,7 @@ AliFMDBaseDigitizer::AliFMDBaseDigitizer()
          AliFMDMap::kMaxRings, 
          AliFMDMap::kMaxSectors, 
          AliFMDMap::kMaxStrips),
-    fShapingTime(0)
+    fShapingTime(6)
 {
   // Default ctor - don't use it
 }
@@ -240,7 +240,7 @@ AliFMDBaseDigitizer::AliFMDBaseDigitizer(AliRunDigitizer* manager)
          AliFMDMap::kMaxRings, 
          AliFMDMap::kMaxSectors, 
          AliFMDMap::kMaxStrips), 
-    fShapingTime(0)
+    fShapingTime(6)
 {
   // Normal CTOR
   AliFMDDebug(1, (" processed"));
@@ -387,7 +387,7 @@ AliFMDBaseDigitizer::DigitizeHits(AliFMD* fmd) const
   //
   AliFMDGeometry* geometry = AliFMDGeometry::Instance();
   
-  TArrayI counts(3);
+  TArrayI counts(4);
   for (UShort_t detector=1; detector <= 3; detector++) {
     AliFMDDebug(5, ("Processing hits in FMD%d", detector));
     // Get pointer to subdetector 
@@ -422,7 +422,7 @@ AliFMDBaseDigitizer::DigitizeHits(AliFMD* fmd) const
          last = edep;
          AddDigit(fmd, detector, ring, sector, strip, edep, 
                   UShort_t(counts[0]), Short_t(counts[1]), 
-                  Short_t(counts[2]));
+                  Short_t(counts[2]), Short_t(counts[3]));
          AliFMDDebug(10, ("   Adding digit in FMD%d%c[%2d,%3d]=%d", 
                          detector,ring,sector,strip,counts[0]));
 #if 0
@@ -495,7 +495,7 @@ AliFMDBaseDigitizer::ConvertToCount(Float_t   edep,
     maxAdc = 1023;
   }
   UShort_t rate           = param->GetSampleRate(detector,ring,sector,strip);
-  if (rate < 1 || rate > 3) rate = 1;
+  if (rate < 1 || rate > 4) rate = 1;
   
   // In case we don't oversample, just return the end value. 
   if (rate == 1) {
@@ -511,7 +511,7 @@ AliFMDBaseDigitizer::ConvertToCount(Float_t   edep,
   // Create a pedestal 
   Float_t b = fShapingTime;
   for (Ssiz_t i = 0; i < rate;  i++) {
-    Float_t t  = Float_t(i) / rate;
+    Float_t t  = Float_t(i) / rate + 1./rate;
     Float_t s  = edep + (last - edep) * TMath::Exp(-b * t);
     Float_t a  = Int_t(s * convF + ped);
     if (a < 0) a = 0;
index 9ae35641bd9e5a55c24556e1500c5339e67f0f9c..9e904fe5e6adddc081575489d2f0df78e92ad85a 100644 (file)
@@ -230,7 +230,8 @@ protected:
                            Float_t  /* edep     */, 
                            UShort_t /* count1   */, 
                            Short_t  /* count2   */, 
-                           Short_t  /* count3   */) const {}
+                           Short_t  /* count3   */,
+                           Short_t  /* count4   */) const {}
 
   AliRunLoader* fRunLoader;       //! Run loader
   AliFMDEdepMap fEdep;             // Cache of Energy from hits 
index 970c1ffeb0303f2f66b63dd981d95fa6ffbf413e..570abd7d0bf36051dd0f285f92ccfb7c31f517e7 100644 (file)
@@ -63,7 +63,7 @@ AliFMDCalibFaker::AliFMDCalibFaker(Int_t mask, const char* loc)
     fPedestalMin(20),
     fPedestalMax(30), 
     fDeadChance(0),
-    fRate(1),
+    fRate(4),
     fZeroThreshold(0),
     fRunMin(0),
     fRunMax(10),
index 178ca82ef9f1640f39f27c5442d77c4b7a50b09f..2de6a9a039a81f75cc15a9af3c4d207e6ef90770 100644 (file)
 
 //====================================================================
 ClassImp(AliFMDDigit)
+#if 0
+; // Here to make Emacs happy
+#endif
 
 //____________________________________________________________________
 AliFMDDigit::AliFMDDigit()
   : fCount1(0),
     fCount2(-1),
-    fCount3(-1)
+    fCount3(-1),
+    fCount4(-1)
 {
   // CTOR
 }
@@ -85,11 +89,13 @@ AliFMDDigit::AliFMDDigit(UShort_t detector,
                         UShort_t strip, 
                         UShort_t count1,
                         Short_t  count2, 
-                        Short_t  count3)
+                        Short_t  count3,
+                        Short_t  count4)
   : AliFMDBaseDigit(detector, ring, sector, strip), 
     fCount1(count1),
     fCount2(count2),
-    fCount3(count3)
+    fCount3(count3), 
+    fCount4(count4)
 {
   //
   // Creates a real data digit object
@@ -122,8 +128,8 @@ AliFMDDigit::Print(Option_t* /* option*/) const
   // Print digit to standard out 
   AliFMDBaseDigit::Print();
   cout << "\t" 
-       << fCount1 << " (+ " << fCount2 << " + " << fCount2 << ") = " 
-       << Counts() << endl;
+       << fCount1 << " (" << fCount2 << "," << fCount3 << "," << fCount4 
+       << ") = " << Counts() << endl;
 }
 
 //____________________________________________________________________
index edb273a6b439f5c13c6113d130832758b83251ec..b865d59202f2e90f7081300f38c265c32697b7e8 100644 (file)
@@ -40,7 +40,8 @@ public:
              UShort_t strip=0, 
              UShort_t count=0, 
              Short_t  count2=-1, 
-             Short_t  count3=-1);
+             Short_t  count3=-1, 
+             Short_t  count4=-1);
   /** DTOR */
   virtual ~AliFMDDigit() {}
   /** @param i # of sample to get 
@@ -52,6 +53,8 @@ public:
   Short_t  Count2()                const { return fCount2;   }
   /** @return ADC count (third sample, or -1 if not used) */
   Short_t  Count3()                const { return fCount3;   }
+  /** @return ADC count (third sample, or -1 if not used) */
+  Short_t  Count4()                const { return fCount4;   }
   /** @return Canonical ADC counts */
   UShort_t Counts()                const;
   /** Print info 
@@ -63,15 +66,17 @@ protected:
   UShort_t fCount1;     // Digital signal 
   Short_t  fCount2;     // Digital signal (-1 if not used)
   Short_t  fCount3;     // Digital signal (-1 if not used)
-  ClassDef(AliFMDDigit,1)     // Normal FMD digit
+  Short_t  fCount4;     // Digital signal (-1 if not used)
+  ClassDef(AliFMDDigit,2)     // Normal FMD digit
 };
 
 inline UShort_t 
 AliFMDDigit::Counts() const 
 {
-  return fCount1 
-    + (fCount2 >= 0 ? fCount2 : 0)
-    + (fCount3 >= 0 ? fCount3 : 0);
+  if (fCount4 >= 0) return fCount3;
+  if (fCount3 >= 0) return fCount2;
+  if (fCount2 >= 0) return fCount2;
+  return fCount1;
 }
 
 inline Int_t
@@ -81,6 +86,7 @@ AliFMDDigit::Count(UShort_t i) const
   case 0: return fCount1;
   case 1: return fCount2;
   case 2: return fCount3;
+  case 3: return fCount4;
   }
   return -1;
 }
index 48e5d5fa68d77521d637b8bd679460f78ddf6f0d..bd421ec7b8dd93318787cee0ffd3aaecbda4a67a 100644 (file)
@@ -327,10 +327,12 @@ AliFMDDigitizer::AddDigit(AliFMD*  fmd,
                          Float_t  /* edep */, 
                          UShort_t count1, 
                          Short_t  count2, 
-                         Short_t  count3) const
+                         Short_t  count3,
+                         Short_t  count4) const
 {
   // Add a digit
-  fmd->AddDigitByFields(detector, ring, sector, strip, count1, count2, count3);
+  fmd->AddDigitByFields(detector, ring, sector, strip, 
+                       count1, count2, count3, count4);
 }
 
 //____________________________________________________________________
@@ -352,6 +354,7 @@ AliFMDDigitizer::CheckDigit(AliFMDDigit*    digit,
   Int_t             integral = counts[0];
   if (counts[1] >= 0) integral += counts[1];
   if (counts[2] >= 0) integral += counts[2];
+  if (counts[3] >= 0) integral += counts[3];
   integral -= Int_t(mean + 2 * width);
   if (integral < 0) integral = 0;
   
index 22ceced04bb6ca690ec45b7ee0df010f54213ae7..384d725c6f01ffae0b46b15aa56ebd9de93e5437 100644 (file)
@@ -70,7 +70,8 @@ protected:
       @param edep     Energy deposited (not used)
       @param count1   ADC count 1
       @param count2   ADC count 2 (-1 if not used)
-      @param count3   ADC count 3 (-1 if not used) */
+      @param count3   ADC count 3 (-1 if not used) 
+      @param count4   ADC count 4 (-1 if not used) */
   virtual void     AddDigit(AliFMD*  fmd,
                            UShort_t detector, 
                            Char_t   ring,
@@ -79,7 +80,8 @@ protected:
                            Float_t  edep, 
                            UShort_t count1, 
                            Short_t  count2, 
-                           Short_t  count3) const;
+                           Short_t  count3,
+                           Short_t  count4) const;
   /** MAke a pedestal
       @param detector Detector #
       @param ring     Ring ID
index dcf488e0b87afe133422205787a0bf2ad8ba9f30..3cd577420c44aa62a5d3662357aad9748fb930ac 100644 (file)
@@ -154,9 +154,11 @@ AliFMDRawReader::ReadAdcs(TClonesArray* array)
                       det, ring, sec, curStr, i));
       new ((*array)[n]) AliFMDDigit(det, ring, sec, curStr, data[i], 
                                    (rate >= 2 ? data[i+1] : 0),
-                                   (rate >= 3 ? data[i+2] : 0));
+                                   (rate >= 3 ? data[i+2] : 0),
+                                   (rate >= 4 ? data[i+3] : 0));
       if (rate >= 2) i++;
       if (rate >= 3) i++;
+      if (rate >= 4) i++;
     }
   }
   return kTRUE;
index 89b4b83c1ad609c8d8cd4996b9fadbea69e633f6..7502788bbb8427b06b0decb95b206c612e7b0061 100644 (file)
@@ -65,6 +65,7 @@ AliFMDRawStream::ReadChannel(UInt_t& ddl, UInt_t& addr,
       AliFMDDebug(30, ("Last is 0x%x, so reading a new word", last));
       next   = Next();
       if(!next){
+       AliFMDDebug(15, ("Read word # %d (!next)", l));
        addr = GetPrevHWAddress();
        ddl  = GetPrevDDLNumber();
        len  = l+1; // Need to add one - l points to last valid index
index acf44e3f102004b15ce631615a2dda39ac6af5f9..055d507d94deb2c630f5378dee8d25d6c99f963e 100644 (file)
@@ -354,8 +354,6 @@ AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
   // load this to subtract a pedestal that was given in a database or
   // something like that. 
 
-  Int_t             counts = 0;
-  Int_t             adc    = 0;
   AliFMDParameters* param  = AliFMDParameters::Instance();
   Float_t           ped    = param->GetPedestal(digit->Detector(), 
                                                digit->Ring(), 
@@ -367,10 +365,11 @@ AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
                                                     digit->Strip());
   AliFMDDebug(15, ("Subtracting pedestal %f from signal %d", 
                   ped, digit->Counts()));
-  if (digit->Count3() > 0)      adc = digit->Count3();
-  else if (digit->Count2() > 0) adc = digit->Count2();
-  else                          adc = digit->Count1();
-  counts = TMath::Max(Int_t(adc - ped), 0);
+  // if (digit->Count3() > 0)      adc = digit->Count3();
+  // else if (digit->Count2() > 0) adc = digit->Count2();
+  // else                          adc = digit->Count1();
+  Int_t adc    = digit->Counts();
+  Int_t counts = TMath::Max(Int_t(adc - ped), 0);
   if (counts < noise * fNoiseFactor) counts = 0;
   if (counts > 0) AliFMDDebug(15, ("Got a hit strip"));
   if (fDiagStep1) fDiagStep1->Fill(adc, counts);
index b7773b334b04ece936582e69d92232847f5db82a..eb3fa44d3c6251de6e063709c92726002c4a5331 100644 (file)
 
 //====================================================================
 ClassImp(AliFMDSDigit)
-
+#if 0
+; // Here to make Emacs happy
+#endif
 //____________________________________________________________________
 AliFMDSDigit::AliFMDSDigit()
   : fEdep(0), 
     fCount1(0),
     fCount2(-1),
-    fCount3(-1)
+    fCount3(-1), 
+    fCount4(-1)
 {
   // cTOR 
 }
@@ -88,12 +91,14 @@ AliFMDSDigit::AliFMDSDigit(UShort_t detector,
                           Float_t  edep,
                           UShort_t count1,
                           Short_t  count2, 
-                          Short_t  count3)
+                          Short_t  count3, 
+                          Short_t  count4)
   : AliFMDBaseDigit(detector, ring, sector, strip), 
     fEdep(edep),
     fCount1(count1),
     fCount2(count2),
-    fCount3(count3)
+    fCount3(count3),
+    fCount4(count4)
 {
   //
   // Creates a real data digit object
@@ -117,8 +122,8 @@ AliFMDSDigit::Print(Option_t* /* option*/) const
   // Print digit to standard out 
   AliFMDBaseDigit::Print();
   cout << "\t" << fEdep << " -> "
-       << fCount1 << " (+ " << fCount2 << " + " << fCount2 << ") = 
-       << Counts() << endl;
+       << fCount1 << " (" << fCount2 << "," << fCount3 << ",
+       << fCount4 << ") = " << Counts() << endl;
 }
 
 //____________________________________________________________________
index 2beb6c9a40ac0f6febf0bc91ba2508142cf61a0d..87003f429febc111107017839ba1cd2710698bd1 100644 (file)
@@ -41,7 +41,8 @@ public:
               Float_t  edep=0,
               UShort_t count=0, 
               Short_t  count2=-1, 
-              Short_t  count3=-1);
+              Short_t  count3=-1,
+              Short_t  count4=-1);
   /** DTOR */
   virtual ~AliFMDSDigit() {}
   /** @return ADC count (first sample) */
@@ -50,6 +51,8 @@ public:
   Short_t  Count2()                const { return fCount2;   }
   /** @return ADC count (third sample, or -1 if not used) */
   Short_t  Count3()                const { return fCount3;   }
+  /** @return ADC count (third sample, or -1 if not used) */
+  Short_t  Count4()                const { return fCount4;   }
   /** @return Canonical ADC counts */
   UShort_t Counts()                const;
   /** @return Energy deposited */
@@ -62,15 +65,17 @@ protected:
   UShort_t fCount1;     // Digital signal 
   Short_t  fCount2;     // Digital signal (-1 if not used)
   Short_t  fCount3;     // Digital signal (-1 if not used)
-  ClassDef(AliFMDSDigit,1)     // Summable FMD digit
+  Short_t  fCount4;     // Digital signal (-1 if not used)
+  ClassDef(AliFMDSDigit,2)     // Summable FMD digit
 };
   
 inline UShort_t 
 AliFMDSDigit::Counts() const 
 {
-  return fCount1 
-    + (fCount2 >= 0 ? fCount2 : 0)
-    + (fCount3 >= 0 ? fCount3 : 0);
+  if (fCount4 >= 0) return fCount3;
+  if (fCount3 >= 0) return fCount2;
+  if (fCount2 >= 0) return fCount2;
+  return fCount1;
 }
 
 
index ea903e17c498628a90ac5bd278dff6b68e3cef72..510cde488f0561e61dbba576aa4748c579423cf8 100644 (file)
@@ -306,11 +306,12 @@ AliFMDSDigitizer::AddDigit(AliFMD*  fmd,
                           Float_t  edep, 
                           UShort_t count1, 
                           Short_t  count2, 
-                          Short_t  count3) const
+                          Short_t  count3,
+                          Short_t  count4) const
 {
   // Add a summable digit
   fmd->AddSDigitByFields(detector, ring, sector, strip, edep, 
-                        count1, count2, count3); 
+                        count1, count2, count3, count4); 
 }
 
 
index 872c95327e6aa84cfd9f50b6aa5ddeb4cd6fd630..c006e8df99c9fa9a3c1c54e183875ddbc49c99ed 100644 (file)
@@ -71,7 +71,8 @@ protected:
                            Float_t  edep, 
                            UShort_t count1, 
                            Short_t  count2, 
-                           Short_t  count3) const;
+                           Short_t  count3,
+                           Short_t  count4) const;
   ClassDef(AliFMDSDigitizer,0) // Make Summable Digits from Hits
 };
 
index cff327acb7676654b29d07aec3d02b1f9bd422d2..2d6fabb3c3e813eb795c3fe9f7781628d91d9844 100644 (file)
@@ -70,13 +70,13 @@ public:
                           eloss.fN-1, eloss.fArray);
     fElossVsPMQ->SetXTitle("p/(mq^{2})=#beta#gamma/q^{2}");
     fElossVsPMQ->SetYTitle("#Delta E/#Delta x / q^{2} [MeV/cm]");
-    fElossVsPMQ->SumW2();
+    fElossVsPMQ->Sumw2();
     fEloss = new TH1D("eloss", "#Delta E/#Delta x / q^{2}", 
                      eloss.fN-1, eloss.fArray);
     fEloss->SetFillColor(2);
     fEloss->SetFillStyle(3001);
     fEloss->SetXTitle("#Delta E/#Delta x / q^{2} [MeV/cm]");
-    fEloss->SumW2();
+    fEloss->Sumw2();
   }
   //__________________________________________________________________
   Bool_t ProcessHit(AliFMDHit* hit, TParticle* p) 
index 628658ae57b3ae12977e9c8c6976b03ef60a3608..7e4cb905610481b5fe4dd1115bc98476b1065ab7 100644 (file)
@@ -26,14 +26,14 @@ MakeCalibration(const char* base="local://$ALICE_ROOT")
 
   gSystem->Load("libFMDutil.so");
   AliFMDCalibFaker f(AliFMDCalibFaker::kAll, 0);
-  f.SetRunRange(0,0);
+  f.SetRunRange(0,999999);
   f.SetGainSeed(AdcPerMip2Gain(60)); // From astrid test beam 
   f.SetThresholdFactor(3);
   f.SetPedestalRange(80,130); // From ASTRID test-beam
   f.SetDeadChance(0);
   f.SetZeroThreshold(0);
   f.SetStripRange(0, 127);
-  f.SetRate(1);
+  f.SetRate(4);
   f.Exec();
 }
 //____________________________________________________________________
index 13b76fcc2b3448683afa0661ba0f0b76e8a00f3d..81815f17ebab0f91f2e4da3f8bc3dd13a754d890 100644 (file)
@@ -10,11 +10,11 @@ void
 RawTest() 
 {
   gRandom->SetSeed(12345);
-  Int_t   sampleRate   = 3;
+  Int_t   sampleRate   = 4;
   Int_t   channelWidth = 128;
   Float_t shapingTime  = 5;
   UInt_t  maxAdc       = (1 << 10);
-  UInt_t  threshold    = (1 << 8);
+  UInt_t  threshold    = 0; // (1 << 8);
   TArrayI outData(sampleRate * channelWidth); 
   
   Float_t lastTotalCharge = 0;
@@ -23,10 +23,10 @@ RawTest()
     Float_t totalCharge = gRandom->Uniform(0, 1);
 
     for (Int_t sample = 0; sample < sampleRate; sample++) {
-      Float_t time   = Float_t(sample) / sampleRate;
+      Float_t time   = Float_t(sample) / sampleRate + 1./sampleRate;
       Float_t charge = (totalCharge + (lastTotalCharge - totalCharge)
                        * TMath::Exp(-shapingTime * time));
-      UInt_t  adc    = UInt_t(maxAdc * charge);      
+      UInt_t  adc    = channel; // UInt_t(maxAdc * charge);      
       outData[channel * sampleRate + sample] = adc;
       if (adc > threshold) ok++;
     }
@@ -37,7 +37,7 @@ RawTest()
            << ")" << std::endl;
   
   { 
-    AliAltroBuffer buffer("FMD_4096.ddl", 1);
+    AliAltroBuffer buffer("FMD_4096.ddl", new AliFMDAltroMapping());
     buffer.WriteDataHeader(kTRUE, kFALSE);
     buffer.WriteChannel(0, 0, 0, outData.fN, outData.fArray, threshold);
     buffer.Flush();
@@ -49,8 +49,8 @@ RawTest()
     std::cerr << "Failed to make AliRawReader" << endl;
     return 0;
   }
-  AliFMDRawStream input(reader, sampleRate);
-  reader->Select(AliFMDParameters::kBaseDDL >> 8);
+  AliFMDRawStream input(reader); // , sampleRate);
+  reader->Select(12); // AliFMDParameters::kBaseDDL >> 8);
   
   Int_t    oldDDL      = -1;
   Int_t    count       = 0;
@@ -63,9 +63,43 @@ RawTest()
   TArrayI counts(10);
   counts.Reset(-1);
   Int_t offset = 0;
+  UInt_t ddl    = 0;
+  UInt_t rate   = 0;
+  UInt_t last   = 0;
+  UInt_t hwaddr = 0;
+  UShort_t data[2048];
 
+  AliFMDParameters* pars = AliFMDParameters::Instance();
   TArrayI inputData(sampleRate * channelWidth); 
+  Bool_t isGood = true;
   while (next) {
+    isGood = input.ReadChannel(ddl, hwaddr, last, data);
+    ddl = 0;
+    
+    std::cout << Form("Read channel %p 0x%x of size %d", ddl, hwaddr, last)
+             << std::endl;
+    UShort_t det, sec, str;
+    Char_t   ring;
+    if (!pars->Hardware2Detector(ddl, hwaddr, det, ring, sec, str)) {
+      std::cerr << Form("Failed to get detector id from DDL %d "
+                       "and hardware address 0x%x", ddl, hwaddr) << std::endl;
+      continue;
+    }
+    rate           = pars->GetSampleRate(det, ring, sec, str);
+    Int_t stripMin = pars->GetMinStrip(det, ring, sec, str);
+    Int_t stripMax = pars->GetMaxStrip(det, ring, sec, str);
+    std::cout << Form("DDL 0x%04x, address 0x%03x maps to FMD%d%c[%2d,%3d]", 
+                     ddl, hwaddr, det, ring, sec, str) << std::endl;
+    
+    // Loop over the `timebins', and make the digits
+    for (size_t i = 0; i < last; i++) {
+      Int_t in  = data[i];
+      Int_t out = outData[channel * sampleRate + sample];
+      std::cout << "[\t" << channel << ",\t" << sample  << "]\t" 
+               << out << "\t" << in << std::flush; 
+      if (out >= threshold && in != out) std::cout << "\tBad" << std::flush;
+    }
+#if 0
     next = input.Next();
 
     if (!next) break;
@@ -81,10 +115,11 @@ RawTest()
              << out << "\t" << in << std::flush; 
     if (out >= threshold && in != out) std::cout << "\tBad" << std::flush;
     std::cout << std::endl;    
+#endif
   }
 
   std::cout << "Read " << count << " values" << std::endl;
-#if 1
+#if 0
   for (Int_t channel = channelWidth - 1; channel > 0; channel--) {
     for (Int_t sample =  sampleRate - 1; sample > 0; sample--) {
       Int_t in  = inputData[channel * sampleRate + sample];
index f03d682023541cde514b437fc69dca8e1a4e4d15..612098ce85a2507155d7ef2ff49b147af6692e82 100644 (file)
@@ -10,6 +10,7 @@ void
 ReadRaw()
 {
   AliCDBManager* cdb = AliCDBManager::Instance();
+  cdb->SetRun(0);
   cdb->SetDefaultStorage("local://$ALICE_ROOT");
   AliLog::SetModuleDebugLevel("FMD", 10);
   AliFMDParameters::Instance()->Init();
diff --git a/FMD/scripts/TestShaping.C b/FMD/scripts/TestShaping.C
new file mode 100644 (file)
index 0000000..d09eb85
--- /dev/null
@@ -0,0 +1,110 @@
+#include <TMath.h>
+#include <TGraph.h>
+#include <TCanvas.h>
+#include <TArrayI.h>
+#include <TRandom.h>
+#include <TLegend.h>
+#include <iostream>
+#include <TH1.h>
+
+void
+convert(UShort_t rate, Int_t adc, Int_t last, TArrayI& counts)
+{
+  if (rate == 1) { 
+    Float_t a = adc;
+    if (a < 0) a = 0;
+    counts[0] = UShort_t(TMath::Min(a, 1023.F));
+    return;
+  }
+  Float_t b = 6;
+  for (Ssiz_t i = 0; i < rate; i++) { 
+    Float_t t = Float_t(i) / rate + 1./rate;
+    Float_t s = adc + (last - adc) * TMath::Exp(-b * t);
+    Float_t a = s;
+    if (a < 0) a = 0;
+    counts[i] = UShort_t(TMath::Min(a, 1023.F));
+    std::cout << " rate=" << rate << "\tadc=" << adc 
+             << "\tcount[" << i << "]=" << counts[i] 
+             << "\ts=" << adc << " + (" << last << " - " 
+             << adc << ") * exp(-" << b << " * " << t 
+             << ")=" << s << std::endl;
+  }
+  return;
+}
+
+TGraph* 
+makeGraph(const TArrayI& adcs, Int_t rate) 
+{
+  Int_t    last = adcs.fArray[0];
+  TArrayI  counts(4);
+  TGraph*  graph = new TGraph(rate * adcs.fN);
+  graph->SetLineColor(rate);
+  graph->SetMarkerColor(rate);
+  graph->SetMarkerStyle(20+rate);
+  graph->SetLineStyle(rate);
+  graph->SetName(Form("rate%d", rate));
+  graph->SetTitle(Form("Rate %d", rate));
+  for (Int_t i = 0; i < adcs.fN; i++) { 
+    counts.Reset(-1);
+    convert(rate, adcs.fArray[i], last, counts);
+    
+    for (Int_t j = 0; j < rate; j++) { 
+      Int_t    idx = (i * rate + j);
+      Double_t x   = (i + (rate > 1 ? Float_t(j+1) / rate-1 : 0));
+      graph->SetPoint(idx, x, counts[j]);
+    }
+    last = counts[rate - 1];
+  }
+  return graph;
+}
+
+void
+TestShaping(int max=4)
+{
+  TArrayI adcs(10);
+  TGraph* orig = new TGraph(adcs.fN);
+  orig->SetName("Original");
+  orig->SetTitle("Original");
+  orig->SetMarkerStyle(25);
+  orig->SetMarkerColor(1);
+  orig->SetMarkerSize(2);
+  orig->SetLineColor(1);
+  for (Int_t i = 0; i < adcs.fN; i++) { 
+    adcs.fArray[i] = Int_t(gRandom->Uniform(0, 1023));
+    orig->SetPoint(i, i, adcs.fArray[i]);
+  }
+
+  TCanvas* c = new TCanvas("c", "c");
+  c->SetFillColor(0);
+  c->SetTopMargin(.02);
+  c->SetRightMargin(.02);
+  
+  TH1* h = new TH1F("frame","frame", adcs.fN+1, -2, adcs.fN);
+  h->SetMinimum(0);
+  h->SetMaximum(1300);
+  h->SetStats(0);
+  h->Draw("");
+  orig->Draw("pl same");
+
+  TLegend* l = new TLegend(adcs.fN*3./4, 1023, adcs.fN, 1300, "", "");
+  l->SetFillColor(0);
+  l->SetBorderSize(1);
+  l->AddEntry(orig, orig->GetTitle(), "lp");
+  
+  for (int i = 1; i <= max; i++) {
+    TGraph* g = makeGraph(adcs, i);
+    g->Draw("pl same");
+    l->AddEntry(g, g->GetTitle(), "lp");
+  }
+  l->Draw();
+  
+  c->Modified();
+  c->Update();
+  c->cd();
+}
+
+  
+  
+    
+
+