]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
- extend debug output for GTU simulation
authorjklein <jklein@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 14 Jul 2011 20:30:23 +0000 (20:30 +0000)
committerjklein <jklein@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 14 Jul 2011 20:30:23 +0000 (20:30 +0000)
- fix order of tracks (crucial for uniquification)
- move pt calculation to gtuParam

TRD/AliTRDgtuParam.cxx
TRD/AliTRDgtuParam.h
TRD/AliTRDgtuSim.cxx
TRD/AliTRDgtuTMU.cxx
TRD/AliTRDtrackGTU.cxx
TRD/AliTRDtrackGTU.h

index 68f1bb8bf8465f2fefd41c0b75e2b5c59394edb7..7243ff09596cf15b3b0390e094ffd01950e7bf4a 100644 (file)
@@ -43,6 +43,8 @@ Bool_t AliTRDgtuParam::fgUseGTUconst = kTRUE;
 // ----- matching windows -----
       Int_t     AliTRDgtuParam::fgDeltaY     = 19;
       Int_t     AliTRDgtuParam::fgDeltaAlpha = 21;
+// ----- reference layers -----
+      Int_t     AliTRDgtuParam::fgRefLayers[] = { 3, 2, 1 };
 
 // ----- Bin widths (granularity) -----
 const Float_t  AliTRDgtuParam::fgkBinWidthY  = 160e-4;
@@ -1194,16 +1196,11 @@ const Bool_t    AliTRDgtuParam::fgZChannelMap[5][16][6][16] = {
 AliTRDgtuParam::AliTRDgtuParam() :
   fVertexSize(20.0),
   fCurrTrackletMask(0),
-  fRefLayers(0x0),
   fMagField(0.5),
   fGeo(0x0)
 {
   // default ctor
   fGeo = new AliTRDgeometry();
-  fRefLayers = new Int_t[fgkNRefLayers];
-  fRefLayers[0] = 3;
-  fRefLayers[1] = 2;
-  fRefLayers[2] = 1;
   for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
     fAki[iLayer] = 0.;
     fBki[iLayer] = 0.;
@@ -1218,7 +1215,6 @@ AliTRDgtuParam::~AliTRDgtuParam()
   // dtor
 
   delete fGeo;
-  delete [] fRefLayers;
 }
 
 AliTRDgtuParam* AliTRDgtuParam::Instance()
@@ -1251,12 +1247,12 @@ Int_t AliTRDgtuParam::GetZSubchannel(Int_t stack, Int_t layer, Int_t zchannel, I
   return fZSubChannel[stack][zchannel][layer][zpos];
 }
 
-Int_t AliTRDgtuParam::GetRefLayer(Int_t refLayerIdx) const
+Int_t AliTRDgtuParam::GetRefLayer(Int_t refLayerIdx)
 {
   // returns the reference layer indexed by refLayerIdx
 
   if (refLayerIdx >= 0 && refLayerIdx < fgkNRefLayers)
-    return fRefLayers[refLayerIdx];
+    return fgRefLayers[refLayerIdx];
   else
     return -1;
 }
@@ -1461,7 +1457,7 @@ Float_t AliTRDgtuParam::GetAki(Int_t k, Int_t i)
   if (fCurrTrackletMask != k)
     GenerateRecoCoefficients(k);
 
-  return fAki[i];
+  return -fAki[i];
 }
 
 Float_t AliTRDgtuParam::GetBki(Int_t k, Int_t i)
@@ -1558,7 +1554,7 @@ Bool_t AliTRDgtuParam::GetIntersectionPoints(Int_t k, Float_t &x1, Float_t &x2)
     return kFALSE;
 }
 
-Int_t AliTRDgtuParam::GetPt(Int_t layerMask, Int_t a, Float_t /* b */, Float_t x1, Float_t x2) const
+Int_t AliTRDgtuParam::GetPt(Int_t layerMask, Int_t a, Float_t /* b */, Float_t x1, Float_t x2, Float_t magField)
 {
   // returns 0.3 * B * 1/a (1/128 GeV/c)
   // a : offset, b : slope (not used)
@@ -1582,7 +1578,7 @@ Int_t AliTRDgtuParam::GetPt(Int_t layerMask, Int_t a, Float_t /* b */, Float_t x
     Int_t layerMaskId = maskIdLut[layerMask];
     Int_t c1 = c1Lut[layerMaskId];
     Int_t c1Ext = c1 << 8;
-    Int_t ptRawStage4 = c1Ext / a;
+    Int_t ptRawStage4 = c1Ext / (a >> 2);
     Int_t ptRawComb4 = ptRawStage4;
     Int_t ptExtComb4 = (ptRawComb4 > 0) ? ptRawComb4 + 33 : ptRawComb4 - 30;
 
@@ -1593,7 +1589,7 @@ Int_t AliTRDgtuParam::GetPt(Int_t layerMask, Int_t a, Float_t /* b */, Float_t x
     Float_t c1 = x1 * x2 / 2. / 10000.; // conversion cm to m
     Float_t r = 0;
     if ( (a >> 1) != 0)
-      r = (0.3 * fMagField / 2. / (fgkBinWidthY/100.)) * (((Int_t) c1) << 8) / (a >> 1); //??? why shift of a?
+      r = (0.3 * magField / 2. / (fgkBinWidthY/100.)) * (((Int_t) c1) << 8) / (a >> 1); //??? why shift of a?
 
     Int_t pt = (Int_t) (2 * r);
     if (pt >= 0)
index 47a30e239d359304a518d3a82163cc3da1e92627..f8fa205912baf8e2bf7098be9920a6863fea57ea 100644 (file)
@@ -50,10 +50,10 @@ class AliTRDgtuParam : public TObject {
   Int_t GetDeltaY() const { return fgDeltaY; }
   Int_t GetDeltaAlpha() const { return fgDeltaAlpha; }
   Int_t GetZSubchannel(Int_t stack, Int_t layer, Int_t zchannel, Int_t zpos) const;
-  Int_t GetRefLayer(Int_t refLayerIdx) const;
+  static Int_t GetRefLayer(Int_t refLayerIdx);
 //  Bool_t GetFitParams(TVectorD &rhs, Int_t k); // const
   Bool_t GetIntersectionPoints(Int_t k, Float_t &x1, Float_t &x2); // const
-  Int_t GetPt(Int_t layerMask, Int_t a, Float_t b, Float_t x1, Float_t x2) const;
+  static Int_t GetPt(Int_t layerMask, Int_t a, Float_t b, Float_t x1, Float_t x2, Float_t magField);
 
   Bool_t IsInZChannel(Int_t stack, Int_t layer, Int_t zchannel, Int_t zpos) const;
 
@@ -99,6 +99,8 @@ class AliTRDgtuParam : public TObject {
   static       Int_t fgDeltaY;         // accepted deviation in y_proj, default: 9
   static       Int_t fgDeltaAlpha;      // accepted deviation in alpha, default: 11
 
+  static       Int_t fgRefLayers[3];    // reference layers for track finding
+
   static       Bool_t fgUseGTUconst;    // use constants as in the GTU for the calculations
                                               // instead of geometry derived quantities
   static const Bool_t fgZChannelMap[5][16][6][16]; // z-channel tables as in GTU
@@ -113,8 +115,6 @@ class AliTRDgtuParam : public TObject {
   Float_t fBki[6]; // coefficients used for the fit, calculated for the current tracklet mask
   Float_t fCki[6]; // coefficients used for the fit, calculated for the current tracklet mask
 
-  Int_t *fRefLayers;           //[fgkNRefLayers] reference layers for track finding
-
   Float_t fMagField;            // magnetic field in T
 
   AliTRDgeometry *fGeo;                //! pointer to the TRD geometry
index 90ae6c29dd03c28b9ad148c5cdc12346e78157f9..8e7dd878a99aaecfc21556da6a655f979e35cb56 100644 (file)
@@ -237,8 +237,10 @@ Bool_t AliTRDgtuSim::RunGTU(AliLoader *loader, AliESDEvent *esd)
            }
            iStackPrev = iStack;
            iSecPrev = iSec;
+           AliDebug(1, Form("now in sec %i, stack %i", iSec, iStack));
        }
-       AliDebug(1, Form("adding tracklet: 0x%08x", trkl->GetTrackletWord()));
+       AliDebug(1, Form("adding tracklet: 0x%08x in sec %i stack %i link %i",
+                        trkl->GetTrackletWord(), trkl->GetDetector() / 30, (trkl->GetDetector() % 30) / 6, trkl->GetHCId() % 12));
        if (fTMU) {
          fTMU->AddTracklet(trkl, iLink);
        }
index 661382db81a5b892685b2e54b51a2051472c9ffe..1e2464e7ecd85a0347b95c8fc5bfb2799ec2fa81 100644 (file)
@@ -195,7 +195,6 @@ Bool_t AliTRDgtuTMU::RunTMU(TList *ListOfTracks, AliESDEvent *esd)
   // ----- track finding -----
   AliDebug(1,"--------- Running tracking units ----------");
   for (Int_t zch = 0; zch < fGtuParam->GetNZChannels(); zch++) {
-    AliDebug(2,Form("Track finder for Zchannel: %i", zch));
     if (!RunTrackFinder(zch, ListOfTracks)) {
       AliError(Form("Track Finder in z-channel %i failed", zch));
       return kFALSE;
@@ -359,7 +358,7 @@ Bool_t AliTRDgtuTMU::RunTrackFinder(Int_t zch, TList* /* ListOfTracks */)
 
    for (Int_t refLayerIdx = 0; refLayerIdx < fGtuParam->GetNRefLayers(); refLayerIdx++) {
      Int_t reflayer = fGtuParam->GetRefLayer(refLayerIdx);
-     AliDebug(5,Form("~~~~~ Reflayer: %i", reflayer));
+     AliDebug(5,Form("Tracking for z-channel: %i, reflayer: %i", zch, reflayer));
 
      ready = kFALSE; // ready if all channels done
 
@@ -406,7 +405,6 @@ Bool_t AliTRDgtuTMU::RunTrackFinder(Int_t zch, TList* /* ListOfTracks */)
        if (0 <= ptrB[reflayer] && ptrB[reflayer] < notr[reflayer])
         trkRB = (AliTRDtrackletGTU*) fZChannelTracklets[reflayer][zch].At(ptrB[reflayer]);
 
-       AliDebug(10,Form("ptrRA: %i, ptrRB: %i", ptrA[reflayer], ptrB[reflayer]));
        yPlus     = trkRA->GetYProj() + fGtuParam->GetDeltaY();
        yMinus    = trkRA->GetYProj() - fGtuParam->GetDeltaY();
        alphaPlus  = trkRA->GetAlpha() + fGtuParam->GetDeltaAlpha();
@@ -507,11 +505,11 @@ Bool_t AliTRDgtuTMU::RunTrackFinder(Int_t zch, TList* /* ListOfTracks */)
         }
        } // end of loop over layers
 
-       AliDebug(5,Form("logic calculation finished, Nhits: %i", nHits));
+       AliDebug(5,Form("logic calculation finished, Nhits: %i %s",
+                      nHits, (nHits >= 4) ? "-> track found" : ""));
 
        if (nHits >= 4) {
         // ----- track registration -----
-        AliDebug(1,"***** TMU: Track found *****");
         AliTRDtrackGTU *track = new AliTRDtrackGTU();
         track->SetSector(fSector);
         track->SetStack(fStack);
@@ -584,7 +582,7 @@ Bool_t AliTRDgtuTMU::RunTrackFinder(Int_t zch, TList* /* ListOfTracks */)
        }
 
        ready = kTRUE;
-       for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
+       for (Int_t layer = fGtuParam->GetNLayers()-1; layer >= 0; layer--) {
 
         bDone[layer] = ptrB[layer] < 0 || ptrB[layer] >= notr[layer];
         ready = ready && bDone[layer];
@@ -593,13 +591,19 @@ Bool_t AliTRDgtuTMU::RunTrackFinder(Int_t zch, TList* /* ListOfTracks */)
           AliError(Form("Invalid increment: %i at ptrA: %i, notr: %i", inc[layer], ptrA[layer], notr[layer]));
 
 //      AliInfo(Form("Shifting layer: %i, notr: %i, ptrA: %i, ptrB: %i, inc: %i", layer, notr[layer], ptrA[layer], ptrB[layer], inc[layer]));
-        AliDebug(10,Form(" -- Layer: %i   %2i(%2i)%s%s   %2i(%2i)%s%s   +%i   %s  (no: %i)",
+        AliDebug(10,Form("    Layer: %i   %2i(%2i, %2i, %4i, %3i)%s%s   %2i(%2i, %2i, %4i, %3i)%s%s   +%i   %s  (no: %i)",
                          layer,
                          ptrA[layer],
                          (0 <= ptrA[layer] && ptrA[layer] < notr[layer]) ? ((AliTRDtrackletGTU*) fZChannelTracklets[layer][zch].At(ptrA[layer]))->GetIndex() : -1,
+                         (0 <= ptrA[layer] && ptrA[layer] < notr[layer]) ? ((AliTRDtrackletGTU*) fZChannelTracklets[layer][zch].At(ptrA[layer]))->GetSubChannel(zch) : -1,
+                         (0 <= ptrA[layer] && ptrA[layer] < notr[layer]) ? ((AliTRDtrackletGTU*) fZChannelTracklets[layer][zch].At(ptrA[layer]))->GetYProj() : -1,
+                         (0 <= ptrA[layer] && ptrA[layer] < notr[layer]) ? ((AliTRDtrackletGTU*) fZChannelTracklets[layer][zch].At(ptrA[layer]))->GetAlpha() : -1,
                          bHitA[layer] ? "*" : " ", bAlignedA[layer] ? "+" : " ",
                          ptrB[layer],
                          (0 <= ptrB[layer] && ptrB[layer] < notr[layer]) ? ((AliTRDtrackletGTU*) fZChannelTracklets[layer][zch].At(ptrB[layer]))->GetIndex() : -1,
+                         (0 <= ptrB[layer] && ptrB[layer] < notr[layer]) ? ((AliTRDtrackletGTU*) fZChannelTracklets[layer][zch].At(ptrB[layer]))->GetSubChannel(zch) : -1,
+                         (0 <= ptrB[layer] && ptrB[layer] < notr[layer]) ? ((AliTRDtrackletGTU*) fZChannelTracklets[layer][zch].At(ptrB[layer]))->GetYProj() : -1,
+                         (0 <= ptrB[layer] && ptrB[layer] < notr[layer]) ? ((AliTRDtrackletGTU*) fZChannelTracklets[layer][zch].At(ptrB[layer]))->GetAlpha() : -1,
                          bHitB[layer] ? "*" : " ", bAlignedB[layer] ? "+" : " ",
                          inc[layer], bDone[layer] ? "done" : "    ", notr[layer]));
         ptrA[layer] += inc[layer];
@@ -683,7 +687,7 @@ Bool_t AliTRDgtuTMU::RunTrackMerging(TList* ListOfTracks)
        TIter trackRefMerged(tracksRefMerged[zch]);
        while (AliTRDtrackGTU *trk = (AliTRDtrackGTU*) trackRefMerged())
          AliDebug(2, Form("track ref layer %i : %i %i %i %i %i %i, y=%i, z_idx=%i",
-                          trk->GetRefLayerIdx(),
+                          AliTRDgtuParam::GetRefLayer(trk->GetRefLayerIdx()),
                           trk->GetTrackletIndex(5),
                           trk->GetTrackletIndex(4),
                           trk->GetTrackletIndex(3),
@@ -696,7 +700,7 @@ Bool_t AliTRDgtuTMU::RunTrackMerging(TList* ListOfTracks)
        TIter trackRefUnique(tracksRefUnique[zch]);
        while (AliTRDtrackGTU *trk = (AliTRDtrackGTU*) trackRefUnique())
          AliDebug(2, Form("track ref layer %i : %i %i %i %i %i %i, y=%i, z_idx=%i",
-                          trk->GetRefLayerIdx(),
+                          AliTRDgtuParam::GetRefLayer(trk->GetRefLayerIdx()),
                           trk->GetTrackletIndex(5),
                           trk->GetTrackletIndex(4),
                           trk->GetTrackletIndex(3),
@@ -712,7 +716,7 @@ Bool_t AliTRDgtuTMU::RunTrackMerging(TList* ListOfTracks)
     do {
        done = kTRUE;
        trkStage0 = 0x0;
-       for (Int_t zch = fGtuParam->GetNZChannels() - 1; zch > -1; zch--) {
+       for (Int_t zch = 0; zch < fGtuParam->GetNZChannels(); zch++) {
            AliTRDtrackGTU *trk = (AliTRDtrackGTU*) tracksRefUnique[zch]->First();
            if (trk == 0) {
                continue;
@@ -722,7 +726,8 @@ Bool_t AliTRDgtuTMU::RunTrackMerging(TList* ListOfTracks)
                minIdx = zch;
                done = kFALSE;
            }
-           else if ( ((trk->GetZChannel() + 3 * trk->GetZSubChannel()) / 2 - 1) < ((trkStage0->GetZChannel() + 3 * trkStage0->GetZSubChannel()) / 2 -1 ) ) {
+           else if ( ((trk->GetZChannel() + 3 * trk->GetZSubChannel()) / 2 - 1) <  ((trkStage0->GetZChannel() + 3 * trkStage0->GetZSubChannel()) / 2 -1 ) ||
+                     (((trk->GetZChannel() + 3 * trk->GetZSubChannel()) / 2 - 1) == ((trkStage0->GetZChannel() + 3 * trkStage0->GetZSubChannel()) / 2 -1 ) && (trk->GetYapprox() < trkStage0->GetYapprox()) ) ) {
                minIdx = zch;
                trkStage0 = trk;
                done = kFALSE;
@@ -737,6 +742,36 @@ Bool_t AliTRDgtuTMU::RunTrackMerging(TList* ListOfTracks)
 
     Uniquifier(tracksZMergedStage0, tracksZUniqueStage0);
 
+    AliDebug(2, "stage 0:");
+    TIter trackZMergedStage0(tracksZMergedStage0);
+    while (AliTRDtrackGTU *trk = (AliTRDtrackGTU*) trackZMergedStage0())
+      AliDebug(2, Form("track ref layer %i : %i %i %i %i %i %i, y=%i, zch=%i idx=%i",
+                      AliTRDgtuParam::GetRefLayer(trk->GetRefLayerIdx()),
+                      trk->GetTrackletIndex(5),
+                      trk->GetTrackletIndex(4),
+                      trk->GetTrackletIndex(3),
+                      trk->GetTrackletIndex(2),
+                      trk->GetTrackletIndex(1),
+                      trk->GetTrackletIndex(0),
+                      trk->GetYapprox() >> 3,
+                      trk->GetZChannel(),
+                      trk->GetZSubChannel()));
+
+    AliDebug(2, "uniquified:");
+    TIter trackZUniqueStage0(tracksZUniqueStage0);
+    while (AliTRDtrackGTU *trk = (AliTRDtrackGTU*) trackZUniqueStage0())
+      AliDebug(2, Form("track ref layer %i : %i %i %i %i %i %i, y=%i, zch=%i idx=%i",
+                      AliTRDgtuParam::GetRefLayer(trk->GetRefLayerIdx()),
+                      trk->GetTrackletIndex(5),
+                      trk->GetTrackletIndex(4),
+                      trk->GetTrackletIndex(3),
+                      trk->GetTrackletIndex(2),
+                      trk->GetTrackletIndex(1),
+                      trk->GetTrackletIndex(0),
+                      trk->GetYapprox() >> 3,
+                      trk->GetZChannel(),
+                      trk->GetZSubChannel()));
+
 // ----- Splitting in z -----
 
     TIter next(tracksZUniqueStage0);
@@ -744,6 +779,23 @@ Bool_t AliTRDgtuTMU::RunTrackMerging(TList* ListOfTracks)
        tracksZSplitted[(trk->GetZChannel() + 3 * (trk->GetZSubChannel() - 1)) % 2]->Add(trk);
     }
 
+    for (Int_t i = 0; i < 2; i++) {
+      AliDebug(2, Form("split %i:", i));
+      TIter trackZSplit(tracksZSplitted[i]);
+      while (AliTRDtrackGTU *trk = (AliTRDtrackGTU*) trackZSplit())
+       AliDebug(2, Form("track ref layer %i : %i %i %i %i %i %i, y=%i, zch=%i idx=%i",
+                        AliTRDgtuParam::GetRefLayer(trk->GetRefLayerIdx()),
+                        trk->GetTrackletIndex(5),
+                        trk->GetTrackletIndex(4),
+                        trk->GetTrackletIndex(3),
+                        trk->GetTrackletIndex(2),
+                        trk->GetTrackletIndex(1),
+                        trk->GetTrackletIndex(0),
+                        trk->GetYapprox() >> 3,
+                        trk->GetZChannel(),
+                        trk->GetZSubChannel()));
+    }
+
 // ----- Merging in zchanels - 2nd stage -----
 
     do {
@@ -759,7 +811,8 @@ Bool_t AliTRDgtuTMU::RunTrackMerging(TList* ListOfTracks)
                minIdx = i;
                done = kFALSE;
            }
-           else if ( ((trk->GetZChannel() + 3 * (trk->GetZSubChannel() - 1)) / 2) < ((trkStage0->GetZChannel() + 3 * (trkStage0->GetZSubChannel() - 1)) / 2) ) {
+           else if ( (((trk->GetZChannel() + 3 * (trk->GetZSubChannel() - 1)) / 2) <  ((trkStage0->GetZChannel() + 3 * (trkStage0->GetZSubChannel() - 1)) / 2)) ||
+                     ((((trk->GetZChannel() + 3 * (trk->GetZSubChannel() - 1)) / 2) == ((trkStage0->GetZChannel() + 3 * (trkStage0->GetZSubChannel() - 1)) / 2)) && (trk->GetYapprox() < trkStage0->GetYapprox()) ) ) {
                minIdx = i;
                trkStage0 = trk;
                done = kFALSE;
@@ -774,6 +827,36 @@ Bool_t AliTRDgtuTMU::RunTrackMerging(TList* ListOfTracks)
 
     Uniquifier(tracksZMergedStage1, ListOfTracks);
 
+    AliDebug(2, "merged:");
+    TIter trackZMergedStage1(tracksZMergedStage1);
+    while (AliTRDtrackGTU *trk = (AliTRDtrackGTU*) trackZMergedStage1())
+      AliDebug(2, Form("track ref layer %i : %i %i %i %i %i %i, y=%i, zch=%i idx=%i",
+                      AliTRDgtuParam::GetRefLayer(trk->GetRefLayerIdx()),
+                      trk->GetTrackletIndex(5),
+                      trk->GetTrackletIndex(4),
+                      trk->GetTrackletIndex(3),
+                      trk->GetTrackletIndex(2),
+                      trk->GetTrackletIndex(1),
+                      trk->GetTrackletIndex(0),
+                      trk->GetYapprox() >> 3,
+                      trk->GetZChannel(),
+                      trk->GetZSubChannel()));
+
+    AliDebug(2, "uniquified:");
+    TIter track(ListOfTracks);
+    while (AliTRDtrackGTU *trk = (AliTRDtrackGTU*) track())
+      AliDebug(2, Form("track ref layer %i : %i %i %i %i %i %i, y=%i, zch=%i idx=%i",
+                      AliTRDgtuParam::GetRefLayer(trk->GetRefLayerIdx()),
+                      trk->GetTrackletIndex(5),
+                      trk->GetTrackletIndex(4),
+                      trk->GetTrackletIndex(3),
+                      trk->GetTrackletIndex(2),
+                      trk->GetTrackletIndex(1),
+                      trk->GetTrackletIndex(0),
+                      trk->GetYapprox() >> 3,
+                      trk->GetZChannel(),
+                      trk->GetZSubChannel()));
+
     // cleaning up
     for (Int_t zch = 0; zch < fGtuParam->GetNZChannels(); zch++) {
       delete tracksRefMerged[zch];
@@ -905,7 +988,7 @@ Bool_t AliTRDgtuTMU::CalculateTrackParams(AliTRDtrackGTU *track)
       AliError(Form("Could not get tracklet in layer %i\n", layer));
       continue;
     }
-    AliDebug(10,Form("trk yprime: %i", trk->GetYPrime()));
+    AliDebug(10,Form("  trk yprime: %6i, aki: %6i", trk->GetYPrime(), (Int_t) (2048 * fGtuParam->GetAki(track->GetTrackletMask(), layer))));
     a += (((Int_t) (2048 * fGtuParam->GetAki(track->GetTrackletMask(), layer))) * trk->GetYPrime() + 1) >> 8;
     b += fGtuParam->GetBki(track->GetTrackletMask(), layer) * trk->GetYPrime() * fGtuParam->GetBinWidthY();
     c += fGtuParam->GetCki(track->GetTrackletMask(), layer) * trk->GetYPrime() * fGtuParam->GetBinWidthY();
@@ -914,14 +997,12 @@ Bool_t AliTRDgtuTMU::CalculateTrackParams(AliTRDtrackGTU *track)
       a += 3;
   a = a >> 2;
 
-  fGtuParam->GetIntersectionPoints(track->GetTrackletMask(), x1, x2);
-  AliDebug(10,Form("Intersection points: %f, %f", x1, x2));
-  AliDebug(10,Form("Sum: a = %5i, b = %9.2f, c = %9.2f\n", a, b, c));
   track->SetFitParams(a, b, c);
 
-  Int_t pt = fGtuParam->GetPt(track->GetTrackletMask(), a, b, x1, x2);
-  track->SetPtInt(pt);
-  AliDebug(5,Form("Track parameters: a = %i, b = %f, c = %f, x1 = %f, x2 = %f, pt = %f (trkl mask: %i)", a, b, c, x1, x2, track->GetPt(), track->GetTrackletMask()));
+  fGtuParam->GetIntersectionPoints(track->GetTrackletMask(), x1, x2);
+
+  AliDebug(5,Form("  Track parameters: a = %i, b = %f, c = %f, x1 = %f, x2 = %f, pt = %f (trkl mask: %i)", a, b, c, x1, x2, track->GetPt(), track->GetTrackletMask()));
+
   return kTRUE;
 }
 
index 4eaa2b9fb63d3d3da2b1b1a4b48545eb322e37da..4a7a52ab0dda2e682f243629676ebf21d34bdb6e 100644 (file)
@@ -42,7 +42,6 @@ AliTRDtrackGTU::AliTRDtrackGTU() :
   TObject(),
   fStack(-1),
   fSector(-1),
-  fPt(0),
   fPID(0),
   fTracklets(0x0),
   fTrackletMask(0),
@@ -137,14 +136,14 @@ Int_t AliTRDtrackGTU::GetZSubChannel()
 
 Int_t AliTRDtrackGTU::GetYapprox()
 {
-// returns an approximated y-position for the track
+  // returns an approximated y-position for the track
+  // taken from the projected y-position of the tracklet in the reference layer
+  // in which the track was found
 
-  for (Int_t layer = 0; layer < AliTRDgtuParam::GetNLayers(); layer++)
-  {
-    if (IsTrackletInLayer(layer))
-      return ((AliTRDtrackletGTU*) (*fTracklets)[layer])->GetYProj();
-  }
-  return 0;
+  if ((fRefLayerIdx > -1) && (fRefLayerIdx < AliTRDgtuParam::GetNRefLayers()))
+    return ((AliTRDtrackletGTU*) (*fTracklets)[AliTRDgtuParam::GetRefLayer(fRefLayerIdx)])->GetYProj();
+  else
+    return 0;
 }
 
 AliESDTrdTrack* AliTRDtrackGTU::CreateTrdTrack() const
index 812c881ddf224b9c186a05e4413c2cc6d98f05b5..3847646d31b6a00c441e96fd469781d372860a69 100644 (file)
@@ -22,8 +22,8 @@ class AliTRDtrackGTU : public TObject {
   ~AliTRDtrackGTU();
 
 // ----- Track properties
-  Int_t    GetPtInt() const { return fPt; }
-  Float_t  GetPt() const { return (Float_t) fPt / 128.; }
+  Int_t    GetPtInt() const { return AliTRDgtuParam::GetPt(fTrackletMask, (Int_t) this->GetA(), 0, 0, 0, 0); }
+  Float_t  GetPt() const { return (Float_t) this->GetPtInt() / 128.; }
   Int_t    GetPID() const { return fPID; }
   Int_t    GetSector() const { return fSector; }
   Int_t    GetStack() const { return fStack; }
@@ -44,6 +44,7 @@ class AliTRDtrackGTU : public TObject {
   Float_t GetC() const { return fC; }
   Int_t GetZChannel() const { return fZChannel; }
   Int_t GetZSubChannel();
+  Int_t GetRefLayer() const { return AliTRDgtuParam::GetRefLayer(fRefLayerIdx); }
   Int_t GetRefLayerIdx() const { return fRefLayerIdx; }
   Int_t GetYapprox();
 
@@ -52,7 +53,6 @@ class AliTRDtrackGTU : public TObject {
 
   void SetStack(Int_t stack) { fStack = stack; }
   void SetSector(Int_t sector) { fSector = sector; }
-  void SetPtInt(Int_t pt) { fPt = pt; }
   void SetPID(Int_t pid) { fPID = pid; }
 
   void SetZChannel(Int_t zch) { fZChannel = zch; }
@@ -68,7 +68,6 @@ class AliTRDtrackGTU : public TObject {
   Int_t fStack; // TRD stack to which this track belongs
   Int_t fSector; // sector in which the track was found
 
-  Int_t fPt; // pt in integer representation
   Int_t fPID; // PID calculated from tracklet PID
 
   TClonesArray *fTracklets; // array holding the tracklets composing this track