ROOT logo
/**************************************************************************
 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
 *                                                                        *
 * Author: The ALICE Off-line Project.                                    *
 * Contributors are mentioned in the code where appropriate.              *
 *                                                                        *
 * Permission to use, copy, modify and distribute this software and its   *
 * documentation strictly for non-commercial purposes is hereby granted   *
 * without fee, provided that the above copyright notice appears in all   *
 * copies and that both the copyright notice and this permission notice   *
 * appear in the supporting documentation. The authors make no claims     *
 * about the suitability of this software for any purpose. It is          *
 * provided "as is" without express or implied warranty.                  *
 **************************************************************************/
/* $Id$ */
/** @file    AliFMDv1.cxx
    @author  Christian Holm Christensen <cholm@nbi.dk>
    @date    Mon Mar 27 12:48:51 2006
    @brief   Concrete implementation of FMD detector driver - detailed
    version 
    @ingroup FMD_sim
*/
//____________________________________________________________________
//                                                                          
// Forward Multiplicity Detector based on Silicon wafers. This class
// contains the base procedures for the Forward Multiplicity detector
// Detector consists of 3 sub-detectors FMD1, FMD2, and FMD3, each of
// which has 1 or 2 rings of silicon sensors. 
// This class contains the detailed version of the FMD - that is, hits
// are produced during simulation. 
//                                                                           
// See also the class AliFMD for a more detailed explanation of the
// various componets. 
//
#include <TVirtualMC.h>		// ROOT_TVirtualMC
#include <AliRun.h>		// ALIRUN_H
#include <AliMC.h>		// ALIMC_H
// #include <AliLog.h>		// ALILOG_H
#include "AliFMDDebug.h" // Better debug macros
#include "AliFMDv1.h"		// ALIFMDV1_H
// #include "AliFMDGeometryBuilder.h"
#include "AliFMDGeometry.h"
#include "AliFMDDetector.h"
#include "AliFMDRing.h"
#include <TParticlePDG.h>
#include <TDatabasePDG.h>
#include "AliFMDHit.h"

//____________________________________________________________________
ClassImp(AliFMDv1)
#if 0
  ; // This is here to keep Emacs for indenting the next line
#endif


//____________________________________________________________________
Bool_t
AliFMDv1::VMC2FMD(TLorentzVector& v, UShort_t& detector,
		  Char_t& ring, UShort_t& sector, UShort_t& strip) const
{
  // Convert VMC coordinates to detector coordinates 
  TVirtualMC*      mc  = TVirtualMC::GetMC();
  AliFMDGeometry*  fmd = AliFMDGeometry::Instance();

  // Get track position
  mc->TrackPosition(v);
  Int_t moduleno; mc->CurrentVolOffID(fmd->GetModuleOff(), moduleno);
  Int_t iring;    mc->CurrentVolOffID(fmd->GetRingOff(), iring);   
  ring = Char_t(iring);
  Int_t det;      mc->CurrentVolOffID(fmd->GetDetectorOff(), det); 
  detector = det;
  

  // Get the ring geometry
  //Int_t     nsec = fmd->GetDetector(detector)->GetRing(ring)->GetNSectors();
  Int_t     nstr  = fmd->GetDetector(detector)->GetRing(ring)->GetNStrips();
  Double_t  lowr  = fmd->GetDetector(detector)->GetRing(ring)->GetMinR();
  Double_t  theta = fmd->GetDetector(detector)->GetRing(ring)->GetTheta();
  Double_t  pitch = fmd->GetDetector(detector)->GetRing(ring)->GetPitch();

  // Figure out the strip number
  Double_t r     = TMath::Sqrt(v.X() * v.X() + v.Y() * v.Y());
  Int_t    str   = Int_t((r - lowr) / pitch);
  if (str < 0 || str >= nstr) return kFALSE;
  strip          = str;

  // Figure out the sector number
  Double_t phi    = TMath::ATan2(v.Y(), v.X()) * 180. / TMath::Pi();
  if (phi < 0) phi = 360. + phi;
  Double_t t      = phi - 2 * moduleno * theta;
  sector          = 2 * moduleno;
  if (t < 0 || t > 2 * theta) return kFALSE;
  else if (t > theta)         sector += 1;

  AliFMDDebug(40, ("<1> Inside an active FMD volume FMD%d%c[%2d,%3d] %s",
		    detector, ring, sector, strip, mc->CurrentVolPath()));
  return kTRUE;
}

//____________________________________________________________________
Bool_t
AliFMDv1::VMC2FMD(Int_t copy, TLorentzVector& v,
		  UShort_t& detector, Char_t& ring,
		  UShort_t& sector, UShort_t& strip) const
{
  // Convert VMC coordinates to detector coordinates 
  TVirtualMC*     mc  = TVirtualMC::GetMC();
  AliFMDGeometry* fmd = AliFMDGeometry::Instance();

  strip = copy - 1;
  Int_t sectordiv; mc->CurrentVolOffID(fmd->GetSectorOff(), sectordiv);
  if (fmd->GetModuleOff() >= 0) {
    Int_t module;    mc->CurrentVolOffID(fmd->GetModuleOff(), module);
    sector = 2 * module + sectordiv;
  }
  else 
    sector = sectordiv;
  AliFMDDebug(30, ("Getting ring volume with offset %d -> %s", 
		    fmd->GetRingOff(), 
		    mc->CurrentVolOffName(fmd->GetRingOff())));
  Int_t iring;     mc->CurrentVolOffID(fmd->GetRingOff(), iring); 
  ring = Char_t(iring);
  Int_t det;       mc->CurrentVolOffID(fmd->GetDetectorOff(), det); 
  detector = det;

  //Double_t  rz  = fmd->GetDetector(detector)->GetRingZ(ring);
  AliFMDDetector* gdet  = fmd->GetDetector(detector);
  AliFMDRing*     gring = gdet->GetRing(ring);
  if (!gring) {
    AliFatal(Form("Ring %c not found (volume was %s at offset %d in path %s)", 
		  ring, 
		  mc->CurrentVolOffName(fmd->GetRingOff()),
		  fmd->GetRingOff(), 
		  mc->CurrentVolPath()));
  }
  Int_t n = gring->GetNSectors();
#if 0
  if (rz < 0) {
    Int_t s = ((n - sector + n / 2) % n) + 1;
    AliFMDDebug(1, ("Recalculating sector to %d (=%d-%d+%d/2%%%d+1 z=%f)",
		     s, n, sector, n, n, rz));
    sector = s;
  }
#endif
  if (sector < 1 || sector > n) {
    AliWarning(Form("sector # %d out of range (0-%d)", sector-1, n-1));
    return kFALSE;
  }
  sector--;
  // Get track position
  mc->TrackPosition(v);
  AliFMDDebug(40, ("<2> Inside an active FMD volume FMD%d%c[%2d,%3d] %s",
		    detector, ring, sector, strip, mc->CurrentVolPath()));

  return kTRUE;
}


//____________________________________________________________________
Bool_t
AliFMDv1::CheckHit(Int_t trackno, Int_t pdg, Float_t absQ, 
		   const TLorentzVector& p, Float_t edep) const
{
  // Check that a hit is good 
  if (AliLog::GetDebugLevel("FMD", "AliFMD") < 5) return kFALSE;
  TVirtualMC* mc   = TVirtualMC::GetMC();
  Double_t mass    = mc->TrackMass();
  Double_t poverm  = (mass == 0 ? 0 : p.P() / mass);
    
  // This `if' is to debug abnormal energy depositions.  We trigger on
  // p/m approx larger than or equal to a MIP, and a large edep - more 
  // than 1 keV - a MIP is 100 eV. 
  if (!(edep > absQ * absQ && poverm > 1)) return kFALSE;
  
  TArrayI procs;
  mc->StepProcesses(procs);
  TString processes;
  for (Int_t ip = 0; ip < procs.fN; ip++) {
    if (ip != 0) processes.Append(",");
    processes.Append(TMCProcessName[procs.fArray[ip]]);
  }
  TDatabasePDG* pdgDB        = TDatabasePDG::Instance();
  TParticlePDG* particleType = pdgDB->GetParticle(pdg);
  TString pname(particleType ? particleType->GetName() : "???");
  TString what;
  if (mc->IsTrackEntering())    what.Append("entering ");
  if (mc->IsTrackExiting())     what.Append("exiting ");
  if (mc->IsTrackInside())      what.Append("inside ");
  if (mc->IsTrackDisappeared()) what.Append("disappeared ");
  if (mc->IsTrackStop())        what.Append("stopped ");
  if (mc->IsNewTrack())         what.Append("new ");
  if (mc->IsTrackAlive())       what.Append("alive ");
  if (mc->IsTrackOut())         what.Append("out ");
      
  Int_t mother = gAlice->GetMCApp()->GetPrimary(trackno);
  AliFMDDebug(15, ("Track # %5d deposits a lot of energy\n" 
		    "  Volume:    %s\n" 
		    "  Momentum:  (%7.4f,%7.4f,%7.4f)\n"
		    "  PDG:       %d (%s)\n" 
		    "  Edep:      %-14.7f keV (mother %d)\n"
		    "  p/m:       %-7.4f/%-7.4f = %-14.7f\n"
		    "  Processes: %s\n"
		    "  What:      %s\n",
		    trackno, mc->CurrentVolPath(), p.X(), p.Y(), p.Z(),
		    pdg, pname.Data(), edep, mother, p.P(), mass, 
		    poverm, processes.Data(), what.Data()));
  return kTRUE;
}


//____________________________________________________________________
void 
AliFMDv1::StepManager()
{
  // Member function that is executed each time a hit is made in the
  // FMD.  None-charged particles are ignored.   Dead tracks  are
  // ignored. 
  //
  // The procedure is as follows: 
  // 
  //   - IF NOT track is alive THEN RETURN ENDIF
  //   - IF NOT particle is charged THEN RETURN ENDIF
  //   - IF NOT volume name is "STRI" or "STRO" THEN RETURN ENDIF 
  //   - Get strip number (volume copy # minus 1)
  //   - Get phi division number (mother volume copy #)
  //   - Get module number (grand-mother volume copy #)
  //   - section # = 2 * module # + phi division # - 1
  //   - Get ring Id from volume name 
  //   - Get detector # from grand-grand-grand-mother volume name 
  //   - Get pointer to sub-detector object. 
  //   - Get track position 
  //   - IF track is entering volume AND track is inside real shape THEN
  //   -   Reset energy deposited 
  //   -   Get track momentum 
  //   -   Get particle ID # 
  ///  - ENDIF
  //   - IF track is inside volume AND inside real shape THEN 
  ///  -   Update energy deposited 
  //   - ENDIF 
  //   - IF track is inside real shape AND (track is leaving volume,
  //         or it died, or it is stopped  THEN
  //   -   Create a hit 
  //   - ENDIF
  //     
  TVirtualMC* mc = TVirtualMC::GetMC();
  if (!mc->IsTrackAlive()) return;
  Double_t absQ = TMath::Abs(mc->TrackCharge());
  if (absQ <= 0) return;
  
  Int_t copy;
  Int_t vol = mc->CurrentVolID(copy);
  AliFMDGeometry*  fmd = AliFMDGeometry::Instance();
  if (!fmd->IsActive(vol)) {
    AliFMDDebug(50, ("Not an FMD volume %d '%s'",vol,mc->CurrentVolName()));
    return;
  }
  TLorentzVector v;
  UShort_t       detector;
  Char_t         ring;
  UShort_t       sector;
  UShort_t       strip;
  
  if (fmd->IsDetailed()) {
    if (!VMC2FMD(copy, v, detector, ring, sector, strip)) return;
  } else {
    if (!VMC2FMD(v, detector, ring, sector, strip)) return;
  }
  TLorentzVector p;
  mc->TrackMomentum(p);
  Int_t    trackno = gAlice->GetMCApp()->GetCurrentTrackNumber();
  Int_t    pdg     = mc->TrackPid();
  Double_t edep    = mc->Edep() * 1000; // keV
  Bool_t   isBad   = CheckHit(trackno, pdg, absQ, p, edep);
  
  // Check that the track is actually within the active area 
  Bool_t entering = mc->IsTrackEntering();
  Bool_t inside   = mc->IsTrackInside();
  Bool_t out      = (mc->IsTrackExiting()|| mc->IsTrackDisappeared()||
		     mc->IsTrackStop());
  // Reset the energy deposition for this track, and update some of
  // our parameters.
  if (entering) {
    AliFMDDebug(15, ("Track # %8d entering active FMD volume %s: "
		      "Edep=%f (%f,%f,%f)", trackno, mc->CurrentVolPath(),
		      edep, v.X(), v.Y(), v.Z()));
    fCurrentP      = p;
    fCurrentV      = v;    
    fCurrentDeltaE = edep;
    fCurrentPdg    = pdg; // mc->IdFromPDG(pdg);
  }
  // If the track is inside, then update the energy deposition
  if (inside && fCurrentDeltaE >= 0) {
    fCurrentDeltaE += edep;
    AliFMDDebug(15, ("Track # %8d inside active FMD volume %s: Edep=%f, "
		      "Accumulated Edep=%f  (%f,%f,%f)", trackno, 
		      mc->CurrentVolPath(), edep, fCurrentDeltaE, 
		      v.X(), v.Y(), v.Z()));
  }
  // The track exits the volume, or it disappeared in the volume, or
  // the track is stopped because it no longer fulfills the cuts
  // defined, then we create a hit. 
  if (out) {
    if (fCurrentDeltaE >= 0) {
      fCurrentDeltaE += edep;
      AliFMDDebug(15, ("Track # %8d exiting active FMD volume %s: Edep=%g, "
			"Accumulated Edep=%g (%f,%f,%f)", trackno, 
			mc->CurrentVolPath(), edep, fCurrentDeltaE, 
			v.X(), v.Y(), v.Z()));
      TVector3 cur(v.Vect());
      cur -= fCurrentV.Vect();
      Double_t len = cur.Mag();
      AliFMDHit* h = 
	AddHitByFields(trackno, detector, ring, sector, strip,
		       fCurrentV.X(),  fCurrentV.Y(), fCurrentV.Z(),
		       fCurrentP.X(),  fCurrentP.Y(), fCurrentP.Z(), 
		       fCurrentDeltaE, fCurrentPdg,   fCurrentV.T(), 
		       len, mc->IsTrackDisappeared()||mc->IsTrackStop());
      // Add a copy 
      if (isBad && fBad) { 
	new ((*fBad)[fBad->GetEntries()]) AliFMDHit(*h);
      }
      // Check the geometry that we can get back the coordinates. 
#ifdef CHECK_TRANS
      Double_t x, y, z;
      fmd->Detector2XYZ(detector, ring, sector, strip, x, y ,z);
      AliFMDDebug(1, ("Hit at (%f,%f,%f), geometry says (%f,%f,%f)", 
		       fCurrentV.X(), fCurrentV.Y(), fCurrentV.Z(), x, y, z));
#endif
    }
    fCurrentDeltaE = -1;
  }
}
//___________________________________________________________________
//
// EOF
//
 AliFMDv1.cxx:1
 AliFMDv1.cxx:2
 AliFMDv1.cxx:3
 AliFMDv1.cxx:4
 AliFMDv1.cxx:5
 AliFMDv1.cxx:6
 AliFMDv1.cxx:7
 AliFMDv1.cxx:8
 AliFMDv1.cxx:9
 AliFMDv1.cxx:10
 AliFMDv1.cxx:11
 AliFMDv1.cxx:12
 AliFMDv1.cxx:13
 AliFMDv1.cxx:14
 AliFMDv1.cxx:15
 AliFMDv1.cxx:16
 AliFMDv1.cxx:17
 AliFMDv1.cxx:18
 AliFMDv1.cxx:19
 AliFMDv1.cxx:20
 AliFMDv1.cxx:21
 AliFMDv1.cxx:22
 AliFMDv1.cxx:23
 AliFMDv1.cxx:24
 AliFMDv1.cxx:25
 AliFMDv1.cxx:26
 AliFMDv1.cxx:27
 AliFMDv1.cxx:28
 AliFMDv1.cxx:29
 AliFMDv1.cxx:30
 AliFMDv1.cxx:31
 AliFMDv1.cxx:32
 AliFMDv1.cxx:33
 AliFMDv1.cxx:34
 AliFMDv1.cxx:35
 AliFMDv1.cxx:36
 AliFMDv1.cxx:37
 AliFMDv1.cxx:38
 AliFMDv1.cxx:39
 AliFMDv1.cxx:40
 AliFMDv1.cxx:41
 AliFMDv1.cxx:42
 AliFMDv1.cxx:43
 AliFMDv1.cxx:44
 AliFMDv1.cxx:45
 AliFMDv1.cxx:46
 AliFMDv1.cxx:47
 AliFMDv1.cxx:48
 AliFMDv1.cxx:49
 AliFMDv1.cxx:50
 AliFMDv1.cxx:51
 AliFMDv1.cxx:52
 AliFMDv1.cxx:53
 AliFMDv1.cxx:54
 AliFMDv1.cxx:55
 AliFMDv1.cxx:56
 AliFMDv1.cxx:57
 AliFMDv1.cxx:58
 AliFMDv1.cxx:59
 AliFMDv1.cxx:60
 AliFMDv1.cxx:61
 AliFMDv1.cxx:62
 AliFMDv1.cxx:63
 AliFMDv1.cxx:64
 AliFMDv1.cxx:65
 AliFMDv1.cxx:66
 AliFMDv1.cxx:67
 AliFMDv1.cxx:68
 AliFMDv1.cxx:69
 AliFMDv1.cxx:70
 AliFMDv1.cxx:71
 AliFMDv1.cxx:72
 AliFMDv1.cxx:73
 AliFMDv1.cxx:74
 AliFMDv1.cxx:75
 AliFMDv1.cxx:76
 AliFMDv1.cxx:77
 AliFMDv1.cxx:78
 AliFMDv1.cxx:79
 AliFMDv1.cxx:80
 AliFMDv1.cxx:81
 AliFMDv1.cxx:82
 AliFMDv1.cxx:83
 AliFMDv1.cxx:84
 AliFMDv1.cxx:85
 AliFMDv1.cxx:86
 AliFMDv1.cxx:87
 AliFMDv1.cxx:88
 AliFMDv1.cxx:89
 AliFMDv1.cxx:90
 AliFMDv1.cxx:91
 AliFMDv1.cxx:92
 AliFMDv1.cxx:93
 AliFMDv1.cxx:94
 AliFMDv1.cxx:95
 AliFMDv1.cxx:96
 AliFMDv1.cxx:97
 AliFMDv1.cxx:98
 AliFMDv1.cxx:99
 AliFMDv1.cxx:100
 AliFMDv1.cxx:101
 AliFMDv1.cxx:102
 AliFMDv1.cxx:103
 AliFMDv1.cxx:104
 AliFMDv1.cxx:105
 AliFMDv1.cxx:106
 AliFMDv1.cxx:107
 AliFMDv1.cxx:108
 AliFMDv1.cxx:109
 AliFMDv1.cxx:110
 AliFMDv1.cxx:111
 AliFMDv1.cxx:112
 AliFMDv1.cxx:113
 AliFMDv1.cxx:114
 AliFMDv1.cxx:115
 AliFMDv1.cxx:116
 AliFMDv1.cxx:117
 AliFMDv1.cxx:118
 AliFMDv1.cxx:119
 AliFMDv1.cxx:120
 AliFMDv1.cxx:121
 AliFMDv1.cxx:122
 AliFMDv1.cxx:123
 AliFMDv1.cxx:124
 AliFMDv1.cxx:125
 AliFMDv1.cxx:126
 AliFMDv1.cxx:127
 AliFMDv1.cxx:128
 AliFMDv1.cxx:129
 AliFMDv1.cxx:130
 AliFMDv1.cxx:131
 AliFMDv1.cxx:132
 AliFMDv1.cxx:133
 AliFMDv1.cxx:134
 AliFMDv1.cxx:135
 AliFMDv1.cxx:136
 AliFMDv1.cxx:137
 AliFMDv1.cxx:138
 AliFMDv1.cxx:139
 AliFMDv1.cxx:140
 AliFMDv1.cxx:141
 AliFMDv1.cxx:142
 AliFMDv1.cxx:143
 AliFMDv1.cxx:144
 AliFMDv1.cxx:145
 AliFMDv1.cxx:146
 AliFMDv1.cxx:147
 AliFMDv1.cxx:148
 AliFMDv1.cxx:149
 AliFMDv1.cxx:150
 AliFMDv1.cxx:151
 AliFMDv1.cxx:152
 AliFMDv1.cxx:153
 AliFMDv1.cxx:154
 AliFMDv1.cxx:155
 AliFMDv1.cxx:156
 AliFMDv1.cxx:157
 AliFMDv1.cxx:158
 AliFMDv1.cxx:159
 AliFMDv1.cxx:160
 AliFMDv1.cxx:161
 AliFMDv1.cxx:162
 AliFMDv1.cxx:163
 AliFMDv1.cxx:164
 AliFMDv1.cxx:165
 AliFMDv1.cxx:166
 AliFMDv1.cxx:167
 AliFMDv1.cxx:168
 AliFMDv1.cxx:169
 AliFMDv1.cxx:170
 AliFMDv1.cxx:171
 AliFMDv1.cxx:172
 AliFMDv1.cxx:173
 AliFMDv1.cxx:174
 AliFMDv1.cxx:175
 AliFMDv1.cxx:176
 AliFMDv1.cxx:177
 AliFMDv1.cxx:178
 AliFMDv1.cxx:179
 AliFMDv1.cxx:180
 AliFMDv1.cxx:181
 AliFMDv1.cxx:182
 AliFMDv1.cxx:183
 AliFMDv1.cxx:184
 AliFMDv1.cxx:185
 AliFMDv1.cxx:186
 AliFMDv1.cxx:187
 AliFMDv1.cxx:188
 AliFMDv1.cxx:189
 AliFMDv1.cxx:190
 AliFMDv1.cxx:191
 AliFMDv1.cxx:192
 AliFMDv1.cxx:193
 AliFMDv1.cxx:194
 AliFMDv1.cxx:195
 AliFMDv1.cxx:196
 AliFMDv1.cxx:197
 AliFMDv1.cxx:198
 AliFMDv1.cxx:199
 AliFMDv1.cxx:200
 AliFMDv1.cxx:201
 AliFMDv1.cxx:202
 AliFMDv1.cxx:203
 AliFMDv1.cxx:204
 AliFMDv1.cxx:205
 AliFMDv1.cxx:206
 AliFMDv1.cxx:207
 AliFMDv1.cxx:208
 AliFMDv1.cxx:209
 AliFMDv1.cxx:210
 AliFMDv1.cxx:211
 AliFMDv1.cxx:212
 AliFMDv1.cxx:213
 AliFMDv1.cxx:214
 AliFMDv1.cxx:215
 AliFMDv1.cxx:216
 AliFMDv1.cxx:217
 AliFMDv1.cxx:218
 AliFMDv1.cxx:219
 AliFMDv1.cxx:220
 AliFMDv1.cxx:221
 AliFMDv1.cxx:222
 AliFMDv1.cxx:223
 AliFMDv1.cxx:224
 AliFMDv1.cxx:225
 AliFMDv1.cxx:226
 AliFMDv1.cxx:227
 AliFMDv1.cxx:228
 AliFMDv1.cxx:229
 AliFMDv1.cxx:230
 AliFMDv1.cxx:231
 AliFMDv1.cxx:232
 AliFMDv1.cxx:233
 AliFMDv1.cxx:234
 AliFMDv1.cxx:235
 AliFMDv1.cxx:236
 AliFMDv1.cxx:237
 AliFMDv1.cxx:238
 AliFMDv1.cxx:239
 AliFMDv1.cxx:240
 AliFMDv1.cxx:241
 AliFMDv1.cxx:242
 AliFMDv1.cxx:243
 AliFMDv1.cxx:244
 AliFMDv1.cxx:245
 AliFMDv1.cxx:246
 AliFMDv1.cxx:247
 AliFMDv1.cxx:248
 AliFMDv1.cxx:249
 AliFMDv1.cxx:250
 AliFMDv1.cxx:251
 AliFMDv1.cxx:252
 AliFMDv1.cxx:253
 AliFMDv1.cxx:254
 AliFMDv1.cxx:255
 AliFMDv1.cxx:256
 AliFMDv1.cxx:257
 AliFMDv1.cxx:258
 AliFMDv1.cxx:259
 AliFMDv1.cxx:260
 AliFMDv1.cxx:261
 AliFMDv1.cxx:262
 AliFMDv1.cxx:263
 AliFMDv1.cxx:264
 AliFMDv1.cxx:265
 AliFMDv1.cxx:266
 AliFMDv1.cxx:267
 AliFMDv1.cxx:268
 AliFMDv1.cxx:269
 AliFMDv1.cxx:270
 AliFMDv1.cxx:271
 AliFMDv1.cxx:272
 AliFMDv1.cxx:273
 AliFMDv1.cxx:274
 AliFMDv1.cxx:275
 AliFMDv1.cxx:276
 AliFMDv1.cxx:277
 AliFMDv1.cxx:278
 AliFMDv1.cxx:279
 AliFMDv1.cxx:280
 AliFMDv1.cxx:281
 AliFMDv1.cxx:282
 AliFMDv1.cxx:283
 AliFMDv1.cxx:284
 AliFMDv1.cxx:285
 AliFMDv1.cxx:286
 AliFMDv1.cxx:287
 AliFMDv1.cxx:288
 AliFMDv1.cxx:289
 AliFMDv1.cxx:290
 AliFMDv1.cxx:291
 AliFMDv1.cxx:292
 AliFMDv1.cxx:293
 AliFMDv1.cxx:294
 AliFMDv1.cxx:295
 AliFMDv1.cxx:296
 AliFMDv1.cxx:297
 AliFMDv1.cxx:298
 AliFMDv1.cxx:299
 AliFMDv1.cxx:300
 AliFMDv1.cxx:301
 AliFMDv1.cxx:302
 AliFMDv1.cxx:303
 AliFMDv1.cxx:304
 AliFMDv1.cxx:305
 AliFMDv1.cxx:306
 AliFMDv1.cxx:307
 AliFMDv1.cxx:308
 AliFMDv1.cxx:309
 AliFMDv1.cxx:310
 AliFMDv1.cxx:311
 AliFMDv1.cxx:312
 AliFMDv1.cxx:313
 AliFMDv1.cxx:314
 AliFMDv1.cxx:315
 AliFMDv1.cxx:316
 AliFMDv1.cxx:317
 AliFMDv1.cxx:318
 AliFMDv1.cxx:319
 AliFMDv1.cxx:320
 AliFMDv1.cxx:321
 AliFMDv1.cxx:322
 AliFMDv1.cxx:323
 AliFMDv1.cxx:324
 AliFMDv1.cxx:325
 AliFMDv1.cxx:326
 AliFMDv1.cxx:327
 AliFMDv1.cxx:328
 AliFMDv1.cxx:329
 AliFMDv1.cxx:330
 AliFMDv1.cxx:331
 AliFMDv1.cxx:332
 AliFMDv1.cxx:333
 AliFMDv1.cxx:334
 AliFMDv1.cxx:335
 AliFMDv1.cxx:336