]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONTrackHitPattern.cxx
Add Config/HighVoltage directory and entry
[u/mrichter/AliRoot.git] / MUON / AliMUONTrackHitPattern.cxx
CommitLineData
7771752e 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
ea94c18b 16/* $Id$ */
17
18
3d1463c8 19//-----------------------------------------------------------------------------
7771752e 20/// \class AliMUONTrackHitPattern
21///
22/// This class propagates tracks to trigger chambers
683cb6c5 23/// searching for matching trigger tracks and fired strips.
7771752e 24///
25/// To each track, a hit pattern for trigger chambers is set.
26/// The hit pattern is a UShort_t with 8 bits used:
fda59e58 27/// <pre>
7771752e 28/// 1 1 0 1 1 1 0 1
29/// | | |
30/// ----------- ------------
31/// chamber: 11 12 13 14 | 11 12 13 14
32/// cathode: bending | non-bending
fda59e58 33/// </pre>
7771752e 34/// The main method is:
683cb6c5 35/// * ExecuteValidation
7771752e 36///
37/// \author Diego Stocco
3d1463c8 38//-----------------------------------------------------------------------------
7771752e 39
40
41#include "AliMUONTrackHitPattern.h"
7ec3b9cf 42
7771752e 43#include "AliMUONConstants.h"
7ec3b9cf 44#include "AliMUONVDigit.h"
45#include "AliMUONDigitMaker.h"
46#include "AliMUONDigitStoreV1.h"
7771752e 47#include "AliMUONGeometryTransformer.h"
7771752e 48#include "AliMUONLocalTrigger.h"
7771752e 49#include "AliMUONLocalTriggerBoard.h"
7ec3b9cf 50#include "AliMUONTrack.h"
51#include "AliMUONTrackExtrap.h"
52#include "AliMUONTrackParam.h"
53#include "AliMUONVTrackStore.h"
54#include "AliMUONVTriggerStore.h"
7771752e 55#include "AliMpPad.h"
7771752e 56#include "AliMpSegmentation.h"
7ec3b9cf 57#include "AliMpVSegmentation.h"
fda59e58 58#include "AliMpDEManager.h"
59#include "AliMUONReconstructor.h"
60#include "AliMUONRecoParam.h"
683cb6c5 61#include "AliMUONTriggerTrack.h"
62#include "AliMUONVTriggerTrackStore.h"
63
64#include "AliMpConstants.h"
ea94c18b 65
66#include "AliMagF.h"
67#include "AliLog.h"
7771752e 68#include "AliTracker.h"
ea94c18b 69
7771752e 70#include <Riostream.h>
7ec3b9cf 71#include <TArrayS.h>
7771752e 72#include <TClonesArray.h>
73#include <TMath.h>
74#include <TMatrixD.h>
683cb6c5 75#include <TROOT.h>
76#include <TDirectory.h>
77#include <TFile.h>
78#include <TSystem.h>
79
80#include <cassert>
7771752e 81
82/// \cond CLASSIMP
83ClassImp(AliMUONTrackHitPattern) // Class implementation in ROOT context
84/// \endcond
85
86
87//______________________________________________________________________________
7ec3b9cf 88AliMUONTrackHitPattern::AliMUONTrackHitPattern(const AliMUONGeometryTransformer& transformer,
89 const AliMUONDigitMaker& digitMaker)
7771752e 90 : TObject(),
7ec3b9cf 91 fTransformer(transformer),
683cb6c5 92 fDigitMaker(digitMaker),
93 fkMaxDistance(99999.)
7771752e 94{
95 /// Default constructor
683cb6c5 96 InitMembers();
7771752e 97
98 // Set magnetic field
99 const AliMagF* kField = AliTracker::GetFieldMap();
100 if (!kField) AliFatal("No field available");
101 AliMUONTrackExtrap::SetField(kField);
7771752e 102}
103
104
105//______________________________________________________________________________
106AliMUONTrackHitPattern::~AliMUONTrackHitPattern(void)
107{
7ec3b9cf 108 /// Destructor
683cb6c5 109 delete fTrigCovariance;
110 for(Int_t cath=0; cath<fgkNcathodes; cath++){
111 delete fPadFired[cath];
112 }
7771752e 113}
114
115
116//______________________________________________________________________________
683cb6c5 117void AliMUONTrackHitPattern::InitMembers()
7771752e 118{
683cb6c5 119 //
120 /// Initialize data members
121 //
122 fDeltaZ = TMath::Abs(AliMUONConstants::DefaultChamberZ(12) - AliMUONConstants::DefaultChamberZ(10));
123
124 const Double_t kTrigNonBendReso = AliMUONConstants::TriggerNonBendingReso();
125 const Double_t kTrigBendReso = AliMUONConstants::TriggerBendingReso();
126 const Double_t kTrigSlopeBendReso = 1.414 * AliMUONConstants::TriggerBendingReso()/fDeltaZ;
127 const Double_t kTrigCovSlopeBend = - kTrigBendReso * kTrigBendReso / fDeltaZ;
128
129 // Covariance matrix 3x3 (X,Y,slopeY) for trigger tracks
130 fTrigCovariance = new TMatrixD(3,3);
131 fTrigCovariance->Zero();
132 (*fTrigCovariance)(0,0) = kTrigNonBendReso * kTrigNonBendReso;
133 (*fTrigCovariance)(1,1) = kTrigBendReso * kTrigBendReso;
134 (*fTrigCovariance)(2,2) = kTrigSlopeBendReso * kTrigSlopeBendReso;
135 (*fTrigCovariance)(1,2) = (*fTrigCovariance)(2,1) = kTrigCovSlopeBend;
136
137 const Int_t kMaxNpads[fgkNcathodes] = {GetMaxX(0)*GetMaxY(0), GetMaxX(1)*GetMaxY(1)};
138 Char_t histoName[40];
139 Char_t *cathCode[fgkNcathodes] = {"bendPlane", "nonBendPlane"};
140
141 const Int_t kNslats = 18;
142
143 for(Int_t cath=0; cath<fgkNcathodes; cath++){
144 sprintf(histoName, "fPadFired%s", cathCode[cath]);
145 fPadFired[cath] = new TH3F(histoName, histoName,
146 fgkNchambers, -0.5, (Float_t)fgkNchambers - 0.5,
147 kNslats, -0.5, (Float_t)kNslats - 0.5,
148 kMaxNpads[cath], -0.5, (Float_t)kMaxNpads[cath] - 0.5);
149 }
150}
fda59e58 151
7ec3b9cf 152
683cb6c5 153//_____________________________________________________________________________
154void AliMUONTrackHitPattern::CheckConstants() const
155{
156/// Check consistence of redefined constants
157
158 assert(fgkNcathodes == AliMpConstants::NofCathodes());
159 assert(fgkNchambers == AliMpConstants::NofTriggerChambers());
160 assert(fgkNplanes == AliMpConstants::NofTriggerChambers() * fgkNcathodes);
161}
162
163
164//______________________________________________________________________________
165void AliMUONTrackHitPattern::ExecuteValidation(AliMUONVTrackStore& trackStore,
166 const AliMUONVTriggerTrackStore& triggerTrackStore,
167 const AliMUONVTriggerStore& triggerStore) const
168{
169 //
170 /// Main method:
171 /// Loops on reco tracks, extrapolates them to trigger chambers
172 /// and searches for matching trigger tracks and digits
173 //
174
175 AliMUONDigitStoreV1 digitStore;
176 TriggerDigits(triggerStore,digitStore);
177
178
179 TIter itTrack(trackStore.CreateIterator());
180 AliMUONTrack* track;
181
182 const Int_t kFirstTrigCh = AliMUONConstants::NTrackingCh();
183
184 while ( ( track = static_cast<AliMUONTrack*>(itTrack()) ) )
185 {
186 AliMUONTrackParam trackParam(*((AliMUONTrackParam*) (track->GetTrackParamAtCluster()->Last())));
187
188 ApplyMCSCorrections(trackParam);
189 AliMUONTrackExtrap::ExtrapToZCov(&trackParam, AliMUONConstants::DefaultChamberZ(kFirstTrigCh)); // extrap to 1st trigger chamber
190
191 AliMUONTriggerTrack *matchedTriggerTrack = MatchTriggerTrack(track, trackParam, triggerTrackStore, triggerStore);
192
193 UShort_t pattern = GetHitPattern(trackParam, matchedTriggerTrack, digitStore);
194 track->SetHitsPatternInTrigCh(pattern);
195 }
196}
197
198
199//______________________________________________________________________________
200AliMUONTriggerTrack *
201AliMUONTrackHitPattern::MatchTriggerTrack(AliMUONTrack* track,
202 AliMUONTrackParam& trackParam,
203 const AliMUONVTriggerTrackStore& triggerTrackStore,
204 const AliMUONVTriggerStore& triggerStore) const
205{
206 //
207 /// Match track with trigger track
208 //
209
210 Int_t matchTrigger = 0;
211 Int_t loTrgNum(-1);
212 Double_t distTriggerTrack[3], sigma2[3];
213 Double_t chi2;
214 Double_t chi2MatchTrigger = 0., minChi2MatchTrigger = 999.;
215 Int_t doubleMatch = -1; // Check if track matches 2 trigger tracks
216 Double_t doubleChi2 = -1.;
217 AliMUONTriggerTrack* doubleTriggerTrack = 0x0;
218 AliMUONTriggerTrack* matchedTriggerTrack = 0x0;
7ec3b9cf 219
683cb6c5 220 const TMatrixD& kParamCov = trackParam.GetCovariances();
7ec3b9cf 221
683cb6c5 222 Double_t xTrack = trackParam.GetNonBendingCoor();
223 Double_t yTrack = trackParam.GetBendingCoor();
224 Double_t ySlopeTrack = trackParam.GetBendingSlope();
225
226 // Covariance matrix 3x3 (X,Y,slopeY) for tracker tracks
227 TMatrixD trackCov(3,3);
228 trackCov.Zero();
229 trackCov(0,0) = kParamCov(0,0);
230 trackCov(1,1) = kParamCov(2,2);
231 trackCov(2,2) = kParamCov(3,3);
232 trackCov(1,2) = kParamCov(2,3);
233 trackCov(2,1) = kParamCov(3,2);
234
235 TMatrixD sumCov(trackCov,TMatrixD::kPlus,*fTrigCovariance);
236
237 Bool_t isCovOK = kTRUE;
238
239 if (sumCov.Determinant() != 0) {
240 sumCov.Invert();
241 } else {
242 AliWarning(" Determinant = 0");
243 isCovOK = kFALSE;
244 sigma2[0] = kParamCov(0,0);
245 sigma2[1] = kParamCov(2,2);
246 sigma2[2] = kParamCov(3,3);
247 // sigma of distributions (trigger-track) X,Y,slopeY
248 const Double_t kDistSigma[3]={AliMUONConstants::TriggerNonBendingReso(),
249 AliMUONConstants::TriggerBendingReso(),
250 1.414 * AliMUONConstants::TriggerBendingReso()/fDeltaZ};
251 for (Int_t iVar = 0; iVar < 3; iVar++) sigma2[iVar] += kDistSigma[iVar] * kDistSigma[iVar];
252 }
fda59e58 253
683cb6c5 254 AliMUONTriggerTrack *triggerTrack;
255 TIter itTriggerTrack(triggerTrackStore.CreateIterator());
256 while ( ( triggerTrack = static_cast<AliMUONTriggerTrack*>(itTriggerTrack() ) ) )
257 {
258 distTriggerTrack[0] = triggerTrack->GetX11() - xTrack;
259 distTriggerTrack[1] = triggerTrack->GetY11() - yTrack;
260 distTriggerTrack[2] = TMath::Tan(triggerTrack->GetThetay()) - ySlopeTrack;
261
262 if(isCovOK){
263 TMatrixD paramDiff(3,1);
264 for(Int_t iVar = 0; iVar < 3; iVar++)
265 paramDiff(iVar,0) = distTriggerTrack[iVar];
266
267 TMatrixD tmp(sumCov,TMatrixD::kMult,paramDiff);
268 TMatrixD chi2M(paramDiff,TMatrixD::kTransposeMult,tmp);
269 chi2 = chi2M(0,0);
270 }
271 else {
272 chi2 = 0.;
273 for (Int_t iVar = 0; iVar < 3; iVar++) chi2 += distTriggerTrack[iVar]*distTriggerTrack[iVar]/sigma2[iVar];
274 }
fda59e58 275
683cb6c5 276 chi2 /= 3.; // Normalized Chi2: 3 degrees of freedom (X,Y,slopeY)
277 if (chi2 < AliMUONReconstructor::GetRecoParam()->GetMaxNormChi2MatchTrigger())
7ec3b9cf 278 {
683cb6c5 279 Bool_t isDoubleTrack = (TMath::Abs(chi2 - minChi2MatchTrigger)<1.);
280 if (chi2 < minChi2MatchTrigger && chi2 < AliMUONReconstructor::GetRecoParam()->GetMaxNormChi2MatchTrigger())
7ec3b9cf 281 {
683cb6c5 282 if(isDoubleTrack)
283 {
284 doubleMatch = loTrgNum;
285 doubleChi2 = chi2MatchTrigger;
286 doubleTriggerTrack = matchedTriggerTrack;
fda59e58 287 }
683cb6c5 288 minChi2MatchTrigger = chi2;
289 chi2MatchTrigger = chi2;
290 loTrgNum = triggerTrack->GetLoTrgNum();
291 matchedTriggerTrack = triggerTrack;
292 AliMUONLocalTrigger* locTrg = triggerStore.FindLocal(loTrgNum);
293 matchTrigger = 1;
294 if(locTrg->LoLpt()>0) matchTrigger = 2;
295 if(locTrg->LoHpt()>0) matchTrigger = 3;
296 }
297 else if(isDoubleTrack)
298 {
299 doubleMatch = triggerTrack->GetLoTrgNum();
300 doubleChi2 = chi2;
7ec3b9cf 301 }
7771752e 302 }
683cb6c5 303 }
304 if(doubleMatch>=0)
305 { // If two trigger tracks match, select the one passing more trigger cuts
306 AliDebug(1, Form("Two candidates found: %i and %i",loTrgNum,doubleMatch));
307 AliMUONLocalTrigger* locTrg1 = triggerStore.FindLocal(doubleMatch);
308 if((locTrg1->LoLpt()>0 && matchTrigger<2) || (locTrg1->LoHpt() && matchTrigger<3))
309 {
310 if(locTrg1->LoHpt()>0) matchTrigger=3;
311 else matchTrigger = 2;
312 loTrgNum = doubleMatch;
313 chi2MatchTrigger = doubleChi2;
314 matchedTriggerTrack = doubleTriggerTrack;
315 }
316 }
317
318 track->SetMatchTrigger(matchTrigger);
319 track->SetLoTrgNum(loTrgNum);
320 track->SetChi2MatchTrigger(chi2MatchTrigger);
321
322 AliMUONLocalTrigger* locTrg = static_cast<AliMUONLocalTrigger*>(triggerStore.FindLocal(loTrgNum));
323
324 if (locTrg)
325 {
326 track->SetLocalTrigger(locTrg->LoCircuit(),
327 locTrg->LoStripX(),
328 locTrg->LoStripY(),
329 locTrg->LoDev(),
330 locTrg->LoLpt(),
331 locTrg->LoHpt());
332 }
333
334 return matchedTriggerTrack;
335}
336
337
338//______________________________________________________________________________
339UShort_t AliMUONTrackHitPattern::GetHitPattern(AliMUONTrackParam &trackParam,
340 AliMUONTriggerTrack* matchedTriggerTrack,
341 AliMUONVDigitStore& digitStore) const
342{
343 //
344 /// Get hit pattern on trigger chambers for the current track
345 //
346 UShort_t pattern = 0;
347 Bool_t isMatch[2];
348 const Int_t kNTrackingCh = AliMUONConstants::NTrackingCh();
349
350 Bool_t patternFromTrigTrack = kFALSE;
351
352
353 // Calculate hit pattern from trigger track
354 if(matchedTriggerTrack){
355 patternFromTrigTrack = PerformTrigTrackMatch(pattern, matchedTriggerTrack, digitStore);
356 }
357
358 if(patternFromTrigTrack) return pattern;
359
360
361 // Calculate hit pattern from tracker track propagation
362 // if hit pattern from trigger track failed
363
364 for(Int_t ch=0; ch<4; ++ch)
365 {
366 Int_t iChamber = kNTrackingCh+ch;
367 AliMUONTrackExtrap::ExtrapToZCov(&trackParam, AliMUONConstants::DefaultChamberZ(iChamber));
368 FindPadMatchingTrack(digitStore, trackParam, isMatch, iChamber);
369 for(Int_t cath=0; cath<2; ++cath)
370 {
371 if(isMatch[cath]) SetBit(pattern, cath, ch);
372 }
373 }
374 return pattern;
375}
376
377
378//______________________________________________________________________________
379void AliMUONTrackHitPattern::SetBit(UShort_t& pattern, Int_t cathode, Int_t chamber) const
380{
381 //
382 /// Set hits pattern
383 //
384 const Int_t kMask[2][4]= {{0x80, 0x40, 0x20, 0x10},
385 {0x08, 0x04, 0x02, 0x01}};
386 pattern |= kMask[cathode][chamber];
387}
388
389
390//______________________________________________________________________________
391void AliMUONTrackHitPattern::AddEffInfo(UShort_t& pattern, Int_t slat, Int_t effType) const
392{
393 //
394 /// Set info on efficiency calculation
395 //
396 pattern += effType << 8;
397 pattern += slat << 10;
398}
399
400
401//______________________________________________________________________________
402void
403AliMUONTrackHitPattern::ApplyMCSCorrections(AliMUONTrackParam& trackParam) const
404{
405 //
406 /// Returns uncertainties on extrapolated position.
407 /// Takes into account Branson plane corrections in the iron wall.
408 //
409
410 const Float_t kZFilterOut = AliMUONConstants::MuonFilterZEnd();
411 const Float_t kFilterThickness = TMath::Abs(kZFilterOut-AliMUONConstants::MuonFilterZBeg()); // cm
412
413 AliMUONTrackExtrap::ExtrapToZCov(&trackParam, kZFilterOut); // Extrap to muon filter end
414 AliMUONTrackExtrap::AddMCSEffect(&trackParam, kFilterThickness, AliMUONConstants::MuonFilterX0()); // Add MCS effects
415 return;
416}
417
418
419//______________________________________________________________________________
420Bool_t
421AliMUONTrackHitPattern::TriggerDigits(const AliMUONVTriggerStore& triggerStore,
422 AliMUONVDigitStore& digitStore) const
423{
424 //
425 /// make (S)Digit for trigger
426 //
427
428 digitStore.Clear();
429
430 AliMUONLocalTrigger* locTrg;
431 TIter next(triggerStore.CreateLocalIterator());
432
433 while ( ( locTrg = static_cast<AliMUONLocalTrigger*>(next()) ) )
434 {
435 if (locTrg->IsNull()) continue;
436
437 TArrayS xyPattern[2];
438 locTrg->GetXPattern(xyPattern[0]);
439 locTrg->GetYPattern(xyPattern[1]);
440
441 // do we need this ? (Ch.F.)
442// for(Int_t cath=0; cath<2; ++cath)
443// {
444// for(Int_t ch=0; ch<4; ++ch)
445// {
446// if(xyPattern[cath][ch]==0) continue;
447// }
448// }
449
450 Int_t nBoard = locTrg->LoCircuit();
451 fDigitMaker.TriggerDigits(nBoard, xyPattern, digitStore);
452 }
453 return kTRUE;
7771752e 454}
455
456
457//______________________________________________________________________________
7ec3b9cf 458void
459AliMUONTrackHitPattern::FindPadMatchingTrack(AliMUONVDigitStore& digitStore,
460 const AliMUONTrackParam& trackParam,
fda59e58 461 Bool_t isMatch[2], Int_t iChamber) const
7771752e 462{
463 //
683cb6c5 464 /// Given the tracker track position, searches for matching digits.
7771752e 465 //
466
467 Float_t minMatchDist[2];
468
7ec3b9cf 469 for(Int_t cath=0; cath<2; ++cath)
470 {
471 isMatch[cath]=kFALSE;
683cb6c5 472 minMatchDist[cath]=fkMaxDistance/10.;
7771752e 473 }
474
7ec3b9cf 475 TIter next(digitStore.CreateIterator());
476 AliMUONVDigit* mDigit;
fda59e58 477
7ec3b9cf 478 while ( ( mDigit = static_cast<AliMUONVDigit*>(next()) ) )
479 {
480 Int_t currDetElemId = mDigit->DetElemId();
fda59e58 481 Int_t currCh = AliMpDEManager::GetChamberId(currDetElemId);
482 if(currCh!=iChamber) continue;
7ec3b9cf 483 Int_t cathode = mDigit->Cathode();
484 Int_t ix = mDigit->PadX();
485 Int_t iy = mDigit->PadY();
486 Float_t xpad, ypad, zpad;
487 const AliMpVSegmentation* seg = AliMpSegmentation::Instance()
488 ->GetMpSegmentation(currDetElemId,AliMp::GetCathodType(cathode));
489
490 AliMpPad pad = seg->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
491 Float_t xlocal1 = pad.Position().X();
492 Float_t ylocal1 = pad.Position().Y();
493 Float_t dpx = pad.Dimensions().X();
494 Float_t dpy = pad.Dimensions().Y();
495 fTransformer.Local2Global(currDetElemId, xlocal1, ylocal1, 0, xpad, ypad, zpad);
496 Float_t matchDist = MinDistanceFromPad(xpad, ypad, zpad, dpx, dpy, trackParam);
497 if(matchDist>minMatchDist[cathode])continue;
498 isMatch[cathode] = kTRUE;
fda59e58 499 if(isMatch[0] && isMatch[1]) break;
7ec3b9cf 500 minMatchDist[cathode] = matchDist;
7771752e 501 }
7771752e 502}
503
504
505//______________________________________________________________________________
7ec3b9cf 506Float_t
507AliMUONTrackHitPattern::MinDistanceFromPad(Float_t xPad, Float_t yPad, Float_t zPad,
508 Float_t dpx, Float_t dpy,
509 const AliMUONTrackParam& trackParam) const
7771752e 510{
511 //
683cb6c5 512 /// Decides if the digit belongs to the tracker track.
7771752e 513 //
7771752e 514
fda59e58 515 AliMUONTrackParam trackParamAtPadZ(trackParam);
516 AliMUONTrackExtrap::ExtrapToZCov(&trackParamAtPadZ, zPad);
517
518 Float_t xTrackAtPad = trackParamAtPadZ.GetNonBendingCoor();
519 Float_t yTrackAtPad = trackParamAtPadZ.GetBendingCoor();
7771752e 520
fda59e58 521 const Float_t kNSigma = AliMUONReconstructor::GetRecoParam()->GetSigmaCutForTrigger();
522
523 const TMatrixD& kCovParam = trackParamAtPadZ.GetCovariances();
524
525 Float_t sigmaX = TMath::Sqrt(kCovParam(0,0));
526 Float_t sigmaY = TMath::Sqrt(kCovParam(2,2));
527
528 Float_t maxDistX = kNSigma * sigmaX; // in cm
529 Float_t maxDistY = kNSigma * sigmaY; // in cm
7771752e 530
531 Float_t deltaX = TMath::Abs(xPad-xTrackAtPad)-dpx;
532 Float_t deltaY = TMath::Abs(yPad-yTrackAtPad)-dpy;
533
683cb6c5 534 Float_t matchDist = fkMaxDistance;
7771752e 535 if(deltaX<=maxDistX && deltaY<=maxDistY) matchDist = TMath::Max(deltaX, deltaY);
fda59e58 536
7771752e 537 return matchDist;
538}
539
540
683cb6c5 541//_____________________________________________________________________________
542Int_t AliMUONTrackHitPattern::FindPadMatchingTrig(AliMUONVDigitStore& digitStore, Int_t &detElemId,
543 Float_t coor[2], Bool_t isMatch[2],
544 TArrayI nboard[2], TArrayF &zRealMatch, Float_t y11) const
545{
546 //
547 /// Check slat and board number of digit matching trigger track
548 //
549
550 enum {kBending, kNonBending};
551
552 Float_t minMatchDist[fgkNcathodes];
553 Int_t padsInCheckArea[fgkNcathodes];
554
555 for(Int_t cath=0; cath<fgkNcathodes; cath++){
556 isMatch[cath] = kFALSE;
557 minMatchDist[cath] = fkMaxDistance/10.;
558 padsInCheckArea[cath] = 0;
559 }
560 Int_t iChamber = AliMpDEManager::GetChamberId(detElemId);
561 Int_t ch = iChamber-10;
562 Float_t oldDeltaZ = AliMUONConstants::DefaultChamberZ(iChamber) - AliMUONConstants::DefaultChamberZ(10);
563 Float_t y = coor[1];
564 Int_t iSlat = detElemId%100;
565 Int_t trigDigitBendPlane = -1;
566 Int_t foundDetElemId = detElemId;
567 Float_t foundZmatch=999.;
568 Float_t yCoorAtPadZ=999.;
569
570 TIter next(digitStore.CreateIterator());
571 AliMUONVDigit* mDigit;
572 Int_t idigit=0;
573
574 while ( ( mDigit = static_cast<AliMUONVDigit*>(next()) ) )
575 {
576 idigit++;
577 Int_t currDetElemId = mDigit->DetElemId();
578 Int_t currCh = AliMpDEManager::GetChamberId(currDetElemId);
579 if(currCh!=iChamber) continue;
580 Int_t currSlat = currDetElemId%100;
581 Int_t slatDiff = TMath::Abs(currSlat-iSlat);
582 if(slatDiff>1 && slatDiff<17) continue; // Check neighbour slats
583 Int_t cathode = mDigit->Cathode();
584 Int_t ix = mDigit->PadX();
585 Int_t iy = mDigit->PadY();
586 Float_t xpad, ypad, zpad;
587 const AliMpVSegmentation* seg = AliMpSegmentation::Instance()
588 ->GetMpSegmentation(currDetElemId,AliMp::GetCathodType(cathode));
589
590 AliMpPad pad = seg->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
591 Float_t xlocal1 = pad.Position().X();
592 Float_t ylocal1 = pad.Position().Y();
593 Float_t dpx = pad.Dimensions().X();
594 Float_t dpy = pad.Dimensions().Y();
595 fTransformer.Local2Global(currDetElemId, xlocal1, ylocal1, 0, xpad, ypad, zpad);
596 AliDebug(2, Form("DetElemId = %i\tCathode = %i\t(x,y) Pad = (%i,%i) = (%.2f,%.2f)\tDim = (%.2f,%.2f)\tTrack = (%.2f,%.2f)\n",currDetElemId,cathode,ix,iy,xpad,ypad,dpx,dpy,coor[0],coor[1]));
597 // searching track intersection with chambers (second approximation)
598 if(ch%2==1){
599 //if(iChamber%2==1){
600 Float_t deltaZ = zpad - zRealMatch[0];
601 y = (coor[1]-y11)*deltaZ/oldDeltaZ + y11;
602 if(TMath::Abs(y-coor[1])>0.1) AliDebug(3, Form("oldDeltaZ = %7.2f newDeltaZ = %7.2f\toldY = %7.2f new y = %7.2f\n",oldDeltaZ,deltaZ,coor[1],y));
603 }
604 Float_t matchDist = PadMatchTrack(xpad, ypad, dpx, dpy, coor[0], y);
605 if(matchDist<fkMaxDistance/2.) padsInCheckArea[cathode]++;
606 if(matchDist>minMatchDist[cathode])continue;
607 isMatch[cathode] = kTRUE;
608 minMatchDist[cathode] = matchDist;
609 foundDetElemId = currDetElemId;
610 foundZmatch=zpad;
611 yCoorAtPadZ=y;
612 if(cathode==kBending) trigDigitBendPlane = idigit;
613 for (Int_t loc=0; loc<pad.GetNofLocations(); loc++){
614 AliMpIntPair location = pad.GetLocation(loc);
615 nboard[cathode][loc] = location.GetFirst();
616 }
617 for(Int_t loc=pad.GetNofLocations(); loc<fgkNlocations; loc++){
618 nboard[cathode][loc]=-1;
619 }
620
621 // Fired pads info
622 Int_t currPair = ix*GetMaxY(cathode) + iy;
623 fPadFired[cathode]->Fill(ch, currSlat, currPair);
624 }
625
626 for(Int_t cath=0; cath<fgkNcathodes; cath++){
627 if(padsInCheckArea[cath]>2) {
628 AliDebug(1, Form("padsInCheckArea[%i] = %i\n",cath,padsInCheckArea[cath]));
629 return -500;
630 }
631 }
632
633 if(isMatch[kBending] || isMatch[kNonBending]){
634 detElemId = foundDetElemId;
635 zRealMatch[ch] = foundZmatch;
636 coor[1] = yCoorAtPadZ;
637 }
638 return trigDigitBendPlane;
639}
640
641//_____________________________________________________________________________
642Float_t AliMUONTrackHitPattern::PadMatchTrack(Float_t xPad, Float_t yPad,
643 Float_t dpx, Float_t dpy,
644 Float_t xTrackAtPad, Float_t yTrackAtPad) const
645{
646 //
647 /// Decides if the digit belongs to the trigger track.
648 //
649
650 Float_t maxDist = 2.;//3. // cm
651 Float_t maxDistCheckArea = 6.; // cm
652
653 Float_t matchDist = fkMaxDistance;
654
655 Float_t deltaX = TMath::Abs(xPad-xTrackAtPad)-dpx;
656 Float_t deltaY = TMath::Abs(yPad-yTrackAtPad)-dpy;
657 Float_t maxDistX = maxDist;
658 Float_t maxDistY = maxDist;
659
660 if(deltaX<=maxDistX && deltaY<=maxDistY) matchDist = TMath::Max(deltaX, deltaY);
661 else if(deltaX<=maxDistCheckArea && deltaY<=maxDistCheckArea) matchDist = fkMaxDistance/5.;
662 return matchDist;
663}
664
665
666//_____________________________________________________________________________
667Int_t AliMUONTrackHitPattern::DetElemIdFromPos(Float_t x, Float_t y,
668 Int_t chamber, Int_t cathode) const
669{
670 //
671 /// Given the (x,y) position in the chamber,
672 /// it returns the corresponding slat
673 //
674
675 Int_t resultingDetElemId = -1;
676 AliMpDEIterator it;
677 Float_t minDist = 999.;
678 for ( it.First(chamber-1); ! it.IsDone(); it.Next() ){
679 Int_t detElemId = it.CurrentDEId();
680 Int_t ich = detElemId/100-10;
681 Float_t tolerance=0.2*((Float_t)ich);
682 Float_t currDist=9999.;
683
684 const AliMpVSegmentation* seg =
685 AliMpSegmentation::Instance()
686 ->GetMpSegmentation(detElemId,AliMp::GetCathodType(cathode));
687 if (!seg) continue;
688
689 Float_t deltax = seg->Dimensions().X();
690 Float_t deltay = seg->Dimensions().Y();
691 Float_t xlocal1 = -deltax;
692 Float_t ylocal1 = -deltay;
693 Float_t xlocal2 = +deltax;
694 Float_t ylocal2 = +deltay;
695 Float_t xg01, yg01, zg1, xg02, yg02, zg2;
696 fTransformer.Local2Global(detElemId, xlocal1, ylocal1, 0, xg01, yg01, zg1);
697 fTransformer.Local2Global(detElemId, xlocal2, ylocal2, 0, xg02, yg02, zg2);
698
699 Float_t xg1 = xg01, xg2 = xg02, yg1 = yg01, yg2 = yg02;
700
701 if(xg01>xg02){
702 xg1 = xg02;
703 xg2 = xg01;
704 }
705 if(yg01>yg02){
706 yg1 = yg02;
707 yg2 = yg01;
708 }
709
710 if(x>=xg1-tolerance && x<=xg2+tolerance && y>=yg1-tolerance && y<=yg2+tolerance){ // takes into account errors in extrapolation
711 if(y<yg1) currDist = yg1-y;
712 else if(y>yg2) currDist = y-yg2;
713 if(currDist<minDist) {
714 resultingDetElemId = detElemId;
715 minDist=currDist;
716 continue;
717 }
718 resultingDetElemId = detElemId;
719 break;
720 }
721 } // loop on detElemId
722 return resultingDetElemId;
723}
724
725
726//_____________________________________________________________________________
727void AliMUONTrackHitPattern::LocalBoardFromPos(Float_t x, Float_t y,
728 Int_t detElemId, Int_t cathode,
729 Int_t localBoard[4]) const
730{
731 //
732 /// Given the (x,y) position in the chamber,
733 /// it returns the corresponding local board
734 //
735
736 for(Int_t loc=0; loc<fgkNlocations; loc++){
737 localBoard[loc]=-1;
738 }
739 Float_t xl, yl, zl;
740 fTransformer.Global2Local(detElemId, x, y, 0, xl, yl, zl);
741 TVector2 pos(xl,yl);
742 const AliMpVSegmentation* seg =
743 AliMpSegmentation::Instance()
744 ->GetMpSegmentation(detElemId,AliMp::GetCathodType(cathode));
745 if (seg){
746 AliMpPad pad = seg->PadByPosition(pos,kFALSE);
747 for (Int_t loc=0; loc<pad.GetNofLocations(); loc++){
748 AliMpIntPair location = pad.GetLocation(loc);
749 localBoard[loc] = location.GetFirst();
750 }
751 }
752}
753
754
755//_____________________________________________________________________________
756Bool_t AliMUONTrackHitPattern::PerformTrigTrackMatch(UShort_t &pattern,
757 AliMUONTriggerTrack *matchedTrigTrack,
758 AliMUONVDigitStore& digitStore) const
7771752e 759{
fda59e58 760 //
683cb6c5 761 /// It searches for matching digits around the trigger track.
fda59e58 762 //
7771752e 763
683cb6c5 764 enum {kBending, kNonBending};
7771752e 765
683cb6c5 766 Int_t chOrder[fgkNchambers] = {0,2,1,3};
767
768 TArrayF zRealMatch(fgkNchambers);
769 TArrayF correctFactor(fgkNcathodes);
770
771 Bool_t isMatch[fgkNcathodes];
772 for(Int_t cath=0; cath<fgkNcathodes; cath++){
773 isMatch[cath] = kFALSE;
774 }
775
776 TArrayF zMeanChamber(fgkNchambers);
777 for(Int_t ch=0; ch<fgkNchambers; ch++){
778 zMeanChamber[ch] = AliMUONConstants::DefaultChamberZ(10+ch);
779 }
780
781 TArrayI digitPerTrack(fgkNcathodes);
782
783 Float_t trackIntersectCh[fgkNchambers][fgkNcathodes];
784
785 TArrayI triggeredDigits;
786 triggeredDigits.Set(fgkNchambers);
787 triggeredDigits.Reset(-1);
788
789 TArrayI trigScheme[fgkNcathodes];
790 TArrayI slatThatTriggered[fgkNcathodes];
791 for(Int_t cath=0; cath<fgkNcathodes; cath++){
792 trigScheme[cath].Set(fgkNchambers);
793 slatThatTriggered[cath].Set(fgkNchambers);
794 }
795
796 Int_t boardThatTriggered[fgkNchambers][fgkNcathodes][fgkNlocations];
797 TArrayI nboard[fgkNcathodes];
798 for(Int_t cath=0; cath<fgkNcathodes; cath++){
799 nboard[cath].Set(fgkNlocations);
800 }
801 Int_t ineffBoard[fgkNlocations];
802 for(Int_t loc=0; loc<fgkNlocations; loc++){
803 ineffBoard[loc] = -1;
804 }
805
806 digitPerTrack.Reset();
807 for(Int_t ch=0; ch<fgkNchambers; ch++){
808 zRealMatch[ch] = zMeanChamber[ch];
809 for(Int_t cath=0; cath<fgkNcathodes; cath++){
810 for(Int_t loc=0; loc<fgkNlocations; loc++){
811 boardThatTriggered[ch][cath][loc]=-1;
812 }
813 }
814 }
815
816 for(Int_t cath=0; cath<fgkNcathodes; cath++){
817 slatThatTriggered[cath].Reset(-1);
818 trigScheme[cath].Reset();
819 }
820
821 Bool_t isClearEvent = kTRUE;
822
823 //Float_t x11 = matchedTrigTrack->GetX11();// x position (info from non-bending plane)
824 Float_t y11 = matchedTrigTrack->GetY11();// y position (info from bending plane)
825 Float_t thetaX = matchedTrigTrack->GetThetax();
826 Float_t thetaY = matchedTrigTrack->GetThetay();
827
828 for(Int_t ch=0; ch<fgkNchambers; ch++) { // chamber loop
829 Int_t currCh = chOrder[ch];
830 AliDebug(2, Form("zMeanChamber[%i] = %.2f\tzRealMatch[0] = %.2f\n",currCh,zMeanChamber[currCh],zRealMatch[0]));
831
832 for(Int_t cath=0; cath<fgkNcathodes; cath++){
833 correctFactor[cath]=1.;
834 }
835 // calculate corrections to trigger track theta
836 if(ch>=1) correctFactor[kNonBending] = zMeanChamber[0]/zRealMatch[0];// corrects x position
837 if(ch>=2) correctFactor[kBending] = (zMeanChamber[2] - zMeanChamber[0]) / (zRealMatch[2] - zRealMatch[0]);// corrects y position
838
839 // searching track intersection with chambers (first approximation)
840 Float_t deltaZ = zMeanChamber[currCh] - zMeanChamber[0];
841 trackIntersectCh[currCh][0] = zMeanChamber[currCh] * TMath::Tan(thetaX) * correctFactor[kNonBending];// x position (info from non-bending plane)
842 trackIntersectCh[currCh][1] = y11 + deltaZ * TMath::Tan(thetaY) * correctFactor[kBending];// y position (info from bending plane)
843 Int_t detElemIdFromTrack = DetElemIdFromPos(trackIntersectCh[currCh][0], trackIntersectCh[currCh][1], 11+currCh, 0);
844 if(detElemIdFromTrack<0) {
845 AliDebug(1, "Warning: trigger track outside trigger chamber\n");
846 continue;
847 }
848
849 triggeredDigits[currCh] = FindPadMatchingTrig(digitStore, detElemIdFromTrack, trackIntersectCh[currCh], isMatch, nboard, zRealMatch, y11);
850
851 // if FindPadMatchingTrig = -500 => too many digits matching pad =>
852 // => Event not clear => Reject track
853 if(triggeredDigits[currCh]<-100){
854 isClearEvent = kFALSE;
855 AliDebug(1, Form("Warning: track = %p (%i) matches many pads. Rejected!\n",(void *)matchedTrigTrack, detElemIdFromTrack));
856 break;
857 }
858
859 for(Int_t cath=0; cath<fgkNcathodes; cath++){
860 if(!isMatch[cath]) continue;
861 SetBit(pattern, cath, currCh);
862 digitPerTrack[cath]++;
863 trigScheme[cath][currCh]++;
864 slatThatTriggered[cath][currCh] = detElemIdFromTrack;
865 for(Int_t loc=0; loc<fgkNlocations; loc++){
866 boardThatTriggered[currCh][cath][loc] = nboard[cath][loc];
867 }
868 }
869 } // end chamber loop
870
871 for(Int_t cath=0; cath<fgkNcathodes; cath++){
872 if(digitPerTrack[cath]<3) isClearEvent = kFALSE;
873 if(!isClearEvent) AliDebug(1, Form("Warning: found %i digits for trigger track cathode %i.\nRejecting event\n", digitPerTrack[cath],cath));
874 }
875
876 if(!isClearEvent) return kFALSE;
877
878 Int_t goodForEff = kBoardEff;
879
880 Int_t ineffSlat = -1;
881 Int_t ineffDetElId = -1;
882 Int_t firstSlat = slatThatTriggered[kBending][0]%100;
883 if(firstSlat<0) firstSlat = slatThatTriggered[kBending][1]%100;
884 Int_t firstBoard = boardThatTriggered[0][kBending][0];
885 if(firstBoard<0) firstBoard = boardThatTriggered[1][kBending][0];
886 for(Int_t ch=0; ch<fgkNchambers; ch++){
887 Bool_t isCurrChIneff = kFALSE;
888 Int_t currSlat = slatThatTriggered[kBending][ch]%100;
889 if(currSlat<0){
890 ineffDetElId = DetElemIdFromPos(trackIntersectCh[ch][0], trackIntersectCh[ch][1], 11+ch, kBending);
891 currSlat = ineffDetElId%100;
892 ineffSlat = currSlat;
893 isCurrChIneff = kTRUE;
894 }
895 if(currSlat!=firstSlat) {
896 AddEffInfo(pattern, 20, kChEff);
897 return kTRUE;
898 }
899 Bool_t atLeastOneLoc=kFALSE;
900 if(isCurrChIneff) LocalBoardFromPos(trackIntersectCh[ch][0], trackIntersectCh[ch][1], ineffDetElId, kBending, ineffBoard);
901 for(Int_t loc=0; loc<fgkNlocations; loc++){
902 Int_t currBoard = boardThatTriggered[ch][kBending][loc];
903 if(isCurrChIneff) currBoard = ineffBoard[loc];
904 if(currBoard==firstBoard){
905 atLeastOneLoc=kTRUE;
906 break;
907 }
908 }
909 if(!atLeastOneLoc) goodForEff = kSlatEff;
910 } // end chamber loop
911
912 AddEffInfo(pattern, firstSlat, goodForEff);
913 return kTRUE;
7771752e 914}
915
916
683cb6c5 917//_____________________________________________________________________________
918void AliMUONTrackHitPattern::UpdateQA() const
7771752e 919{
7ec3b9cf 920 //
683cb6c5 921 /// Save map of fired strips in the QA file
7ec3b9cf 922 //
24754ef8 923
683cb6c5 924 TDirectory *dir = gDirectory;
925
926 TFile *logFile = 0x0;
927 TString logFileName;
928 Bool_t reopenFile = kFALSE;
929 TString baseFileName = "MUON.QA";
930 TString dirName = "EfficiencyRelated";
931
932 TSeqCollection *list = gROOT->GetListOfFiles();
933 Int_t n = list->GetEntries();
934 for(Int_t i=0; i<n; i++) {
935 logFile = (TFile*)list->At(i);
936 logFileName = logFile->GetName();
937 if (logFileName.Contains(baseFileName.Data())) break;
938 logFile = 0x0;
939 }
940
941 if(!logFile) {
942 void * dirp = gSystem->OpenDirectory(gSystem->pwd());
943 const char * name = 0x0;
944 // Add all files matching *pattern* to the chain
945 while((name = gSystem->GetDirEntry(dirp))) {
946 logFileName = name;
947 if (logFileName.Contains(baseFileName.Data())) {
948 logFile = new TFile(logFileName.Data(), "update");
949 AliWarning(Form("%s already stored on disk. Re-opening in update mode.",baseFileName.Data()));
950 break;
951 }
952 logFile = 0x0;
953 }//directory loop
954 reopenFile = kTRUE;
955 }
7ec3b9cf 956
683cb6c5 957 if(logFile){
958 logFile->cd();
959 TDirectory *muonDir = (TDirectory*)logFile->Get("MUON");
960 muonDir->cd();
961 TDirectory *effDir = (TDirectory*)muonDir->Get(dirName.Data());
962 if(!effDir) {
963 effDir = muonDir->mkdir(dirName.Data());
964 }
965 effDir->cd();
966 TH3F *histo = 0x0;
967 for(Int_t cath=0; cath<fgkNcathodes; cath++){
968 histo = (TH3F*) effDir->Get(fPadFired[cath]->GetName());
969 if(!histo) histo = (TH3F*)fPadFired[cath]->Clone();
970 else histo->Add(fPadFired[cath]);
971 histo->Write(histo->GetName(), TObject::kOverwrite);
972 }
973 if(reopenFile){
974 logFile->Close();
975 }
7ec3b9cf 976 }
683cb6c5 977 else AliWarning(Form("Map of strips entering efficiency calculation could not be found in %s", baseFileName.Data()));
978 dir->cd();
7771752e 979}
683cb6c5 980