ROOT logo
// This is an example macro which loops on all ESD events in a file set
// and stores the selected information in form of mini ESD tree.
// It also extracts the corresponding MC information in a parallel MC tree.
//
// Input: galice.root,Kinematics.root,AliESDs.root (read only)
// Output: miniesd.root (recreate)
//
// Usage: faster if compiled
// gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include -I$ALICE_ROOT/TPC");  
// gROOT->LoadMacro("mcminiesd.C+");
// mcminiesd()
//
// Author: Peter Hristov
// e-mail: Peter.Hristov@cern.ch

#if !defined(__CINT__) || defined(__MAKECINT__)

// Root include files
#include <Riostream.h>
#include <TFile.h>
#include <TTree.h>
#include <TBranch.h>
#include <TStopwatch.h>
#include <TObject.h>
#include <TParticle.h>
#include <TLorentzVector.h>
#include <TVector3.h>
#include <TArrayF.h>
#include <TDatabasePDG.h>
#include <TParticlePDG.h>

// AliRoot include files
#include "AliESD.h"
#include "AliESDVertex.h"
#include "AliRunLoader.h"
#include "AliRun.h"
#include "AliStack.h"
#include "AliHeader.h"
#include "AliGenEventHeader.h"
#include "AliTPCtrack.h"
#include "AliTracker.h"

#endif

void selectRecTracks(AliESD* esdOut, Int_t * select, Int_t &nkeep0) {
  // Select also all the reconstructed tracks in case it was not done before
  Int_t nrec = esdOut->GetNumberOfTracks();
  for(Int_t irec=0; irec<nrec; irec++) {
    AliESDtrack * track = esdOut->GetTrack(irec);
    UInt_t label = TMath::Abs(track->GetLabel());
    if (label>=10000000) continue; // Underlying event. 10000000 is the
                                   // value of fkMASKSTEP in AliRunDigitizer
    if (!select[label]) {
      select[label] = 1;
      nkeep0++;
    }
  }
}

void copyGeneralESDInfo(AliESD* esdIn, AliESD* esdOut) {
  // Event and run number
  esdOut->SetEventNumber(esdIn->GetEventNumber());
  esdOut->SetRunNumber(esdIn->GetRunNumber());
  
  // Trigger
  esdOut->SetTriggerMask(esdIn->GetTriggerMask());

  // Magnetic field
  esdOut->SetMagneticField(esdIn->GetMagneticField());

  // Copy ESD vertex
  const AliESDVertex * vtxIn = esdIn->GetVertex();
  AliESDVertex * vtxOut = 0x0;
  if (vtxIn) {
    Double_t pos[3];
    vtxIn->GetXYZ(pos);
    Double_t cov[6];
    vtxIn->GetCovMatrix(cov);
  
    vtxOut = new AliESDVertex(pos,cov,
					     vtxIn->GetChi2(),
					     vtxIn->GetNContributors());
  }
  else
    vtxOut = new AliESDVertex();
  
  esdOut->SetVertex(vtxOut);
}

void selectMiniESD(AliESD* esdIn, AliESD* &esdOut) {
  // This function copies the general ESD information
  // and selects the reconstructed tracks of interest
  // The second argument is a reference to pointer since we 
  // want to operate not with the content, but with the pointer itself

  printf("--------------------\n");
  printf("Selecting data mini ESD\n");
  TStopwatch timer;
  timer.Start();

  // Create the new output ESD
  esdOut = new AliESD();

  // Copy the general information
  copyGeneralESDInfo(esdIn, esdOut);

  // Select tracks
  Int_t ntrk = esdIn->GetNumberOfTracks();
  
  // Loop on tracks
  for (Int_t itrk=0; itrk<ntrk; itrk++) {
    
    AliESDtrack * trackIn = esdIn->GetTrack(itrk);
    UInt_t status=trackIn->GetStatus();

    //select only tracks with TPC or ITS refit
    if ((status&AliESDtrack::kTPCrefit)==0
	&& (status&AliESDtrack::kITSrefit)==0
	) continue;

    //select only tracks with "combined" PID
    if ((status&AliESDtrack::kESDpid)==0) continue;
    
    AliESDtrack * trackOut = new AliESDtrack(*trackIn);
    
    // Reset most of the redundant information
    trackOut->MakeMiniESDtrack();
    
    esdOut->AddTrack(trackOut);
    
    // Do not delete trackIn, it is a second pointer to an
    // object belonging to the TClonesArray
    delete trackOut;
  }

  printf("%d tracks selected\n",esdOut->GetNumberOfTracks());
  esdOut->Print();
  timer.Stop();
  timer.Print();
  printf("--------------------\n");
}

void fixMotherDaughter(TClonesArray& particles, Int_t * map) {
  // Fix mother-daughter relationship
  // using the map of old to new indexes
  Int_t nkeep = particles.GetEntriesFast();
  for (Int_t i=0; i<nkeep; i++) {
    TParticle * part = (TParticle*)particles[i];
    
    // mother
    Int_t mum = part->GetFirstMother();      
    if (mum>-1 && i>map[mum]) 
      part->SetFirstMother(map[mum]);
    
    // old indexes
    Int_t fd = part->GetFirstDaughter();
    Int_t ld = part->GetLastDaughter();
    
    // invalidate daughters
    part->SetFirstDaughter(-1);
    part->SetLastDaughter(-1);
    
    for (Int_t id=TMath::Max(fd,0); id<=ld; id++) {
      if (map[id]>-1) {
	// this daughter was selected
	if(part->GetFirstDaughter() < 0) part->SetFirstDaughter(map[id]);
	part->SetLastDaughter(map[id]);
      }
    }
  }
}

void checkConsistency(TClonesArray& particles) {
  // Check consistency
  // each mother is before the daughters,
  // each daughter from the mother's list 
  // points back to its mother
  Int_t nkeep = particles.GetEntriesFast();
  for (Int_t i=0; i<nkeep; i++) {
    TParticle * part = (TParticle*)particles[i];
    
    Int_t mum = part->GetFirstMother();
    
    if (mum>-1 && i<mum) {
      cout << "Problem: mother " << mum << " after daughter " << i << endl;
    }
    
    Int_t fd = part->GetFirstDaughter();
    Int_t ld = part->GetLastDaughter();
    
    if (fd > ld ) cout << "Problem " << fd << " > " << ld << endl;
    if (fd > -1 && fd < i ) 
      cout << "Problem: mother " << i << " after daughter " << ld << endl;
    
    for (Int_t id=TMath::Max(fd,0); id<=ld; id++) {
      TParticle * daughter = (TParticle*)particles[id];
      if (daughter->GetFirstMother() != i) {
	cout << "Problem "<< i << " not " 
	     << daughter->GetFirstMother()  << endl;
	daughter->Print();
      }
    }
    
  }
}


void kinetree(AliRunLoader* runLoader, AliESD* &esdOut, TClonesArray* &ministack) {

  printf("--------------------\n");
  printf("Selecting MC mini ESD\n");
  TStopwatch timer;
  timer.Start();

  // Particle properties
  TDatabasePDG * pdgdb = TDatabasePDG::Instance();

  // Particle stack
  AliStack * stack = runLoader->Stack();

  Int_t npart = stack->GetNtrack();

  // Particle momentum and vertex. Will be extracted from TParticle
  TLorentzVector momentum;
  TArrayF vertex(3);
  runLoader->GetHeader()->GenEventHeader()->PrimaryVertex(vertex);
  TVector3 dvertex;

  // Counter for selected particles
  Int_t nkeep0 = 0;

  Int_t * select = new Int_t[npart];

  // Loop on particles: physics selection
  for (Int_t ipart=0; ipart<npart; ipart++) {
      
    select[ipart] = 0;

    TParticle * part = stack->Particle(ipart);

    if (!part) {
      cout << "Empty particle " << ipart << endl;
      continue;
    }

    // Particle momentum and vertex of origin
    part->Momentum(momentum);
    dvertex.SetXYZ(part->Vx()-vertex[0],
		   part->Vy()-vertex[1],
		   part->Vz()-vertex[2]);

    // Now select only the "interesting" particles

    // Select all particles from the primary vertex:
    // resonances and products of resonance decays
    if(dvertex.Mag()<0.0001) {
      select[ipart] = 1;
      nkeep0++;
      continue;
    }

    // Reject particles born outside ITS
    if (dvertex.Perp()>100) continue;

    // Reject particles with too low momentum, they don't rich TPC
    if (part->Pt()<0.1) continue;

    // Reject "final state" neutral particles
    // This has to be redone after this loop
    Int_t pdgcode = part->GetPdgCode();
    TParticlePDG * pdgpart = pdgdb->GetParticle(pdgcode);
    Int_t charge  = 0;
    if (pdgpart) charge = (Int_t) pdgpart->Charge();

    if (!charge) {
      if (part->GetFirstDaughter()<0) continue;
    }
      
    // Select the rest of particles except the case
    // particle -> particle + X
    // for example bremstrahlung and delta electrons
    Int_t mumid = part->GetFirstMother();
    TParticle * mother = stack->Particle(mumid);
    if (mother) {
      Int_t mumpdg = mother->GetPdgCode();
      Bool_t skip = kFALSE;
      for (Int_t id=mother->GetFirstDaughter();
	   id<=mother->GetLastDaughter();
	   id++) {
	TParticle * daughter=stack->Particle(id);
	if (!daughter) continue;
	if (mumpdg == daughter->GetPdgCode()) {
	  skip=kTRUE;
	  break;
	}
      }
      if (skip) continue;
    }
    
    // All the criteria are OK, this particle is selected
    select[ipart] = 1;
    nkeep0++;

  } // end loop over particles in the current event

  selectRecTracks(esdOut, select, nkeep0);

  // Container for the selected particles
  TClonesArray * pparticles = new TClonesArray("TParticle",nkeep0);
  TClonesArray &particles = *pparticles;

  // Map of the indexes in the stack and the new ones in the TClonesArray
  Int_t * map = new Int_t[npart];

  // Counter for selected particles
  Int_t nkeep = 0;

  for (Int_t ipart=0; ipart<npart; ipart++) {
      
    map[ipart] = -99; // Reset the map
      
    if (select[ipart]) {

      TParticle * part = stack->Particle(ipart);
      map[ipart] = nkeep;
      new(particles[nkeep++]) TParticle(*part);

    }
  }

  if (nkeep0 != nkeep) printf("Strange %d is not %d\n",nkeep0,nkeep);

  delete [] select;

  // Fix mother-daughter relationship
  fixMotherDaughter(particles,map);

  // Check consistency
  checkConsistency(particles);

  // Now remove the "final state" neutral particles

  TClonesArray * particles1 = new TClonesArray("TParticle",nkeep);
  Int_t * map1 = new Int_t[nkeep];
  Int_t nkeep1 = 0;
    
  // Loop on particles
  for (Int_t ipart=0; ipart<nkeep; ipart++) {
      
    map1[ipart] = -99; // Reset the map

    TParticle * part = (TParticle*)particles[ipart];
    
    // Reject "final state" neutral particles
    // This has to be redone after this loop
    Int_t pdgcode = part->GetPdgCode();
    TParticlePDG * pdgpart = pdgdb->GetParticle(pdgcode);
    Int_t charge  =  0;
    if (pdgpart) charge = (Int_t) pdgpart->Charge();
    
    if (!charge) {
      if (part->GetFirstDaughter()<0) continue;
    }

    // Select the particle
    map1[ipart] = nkeep1;
    TClonesArray &ar = *particles1;
    new(ar[nkeep1++]) TParticle(*part);
  }

  particles.Delete();
  delete pparticles;
  cout << nkeep1 << " particles finally selected" << endl;

  fixMotherDaughter(*particles1,map1);
  checkConsistency(*particles1);

  // Remap the labels of reconstructed tracks

  AliESD * esdNew = new AliESD();

  copyGeneralESDInfo(esdOut,esdNew);

  // Tracks
  Int_t nrec = esdOut->GetNumberOfTracks();
  for(Int_t irec=0; irec<nrec; irec++) {
    AliESDtrack * track = esdOut->GetTrack(irec);
    UInt_t label = TMath::Abs(track->GetLabel());
    if (label<10000000) { // Signal event. 10000000 is the
                          // value of fkMASKSTEP in AliRunDigitizer
      track->SetLabel(TMath::Sign(map1[map[label]],track->GetLabel()));
    }
    esdNew->AddTrack(track);
  }

  delete esdOut;
  esdOut = esdNew;

  ministack = particles1;


  delete [] map;
  delete [] map1;


  timer.Stop();
  timer.Print();
  printf("--------------------\n");
}



AliTPCtrack * particle2track(TParticle * part) {

  // Converts TParticle to AliTPCtrack

  UInt_t index;
  Double_t xx[5];
  Double_t cc[15];
  Double_t xref;
  Double_t alpha;

  // Calculate alpha: the rotation angle of the corresponding TPC sector
  alpha = part->Phi()*180./TMath::Pi();
  if (alpha<0) alpha+= 360.;
  if (alpha>360) alpha -= 360.;

  Int_t sector = (Int_t)(alpha/20.);
  alpha = 10. + 20.*sector;
  alpha /= 180;
  alpha *= TMath::Pi();


  // Reset the covariance matrix
  for (Int_t i=0; i<15; i++) cc[i]=0.;


  // Index
  index = part->GetUniqueID();

  // Get the vertex of origin and the momentum
  TVector3 ver(part->Vx(),part->Vy(),part->Vz());
  TVector3 mom(part->Px(),part->Py(),part->Pz());


  // Rotate to the local (sector) coordinate system
  ver.RotateZ(-alpha);
  mom.RotateZ(-alpha);

  // X of the referense plane
  xref = ver.X();

  // Track parameters
  // fX     = xref         X-coordinate of this track (reference plane)
  // fAlpha = Alpha        Rotation angle the local (TPC sector) 
  // fP0    = YL           Y-coordinate of a track
  // fP1    = ZG           Z-coordinate of a track
  // fP2    = C*x0         x0 is center x in rotated frame
  // fP3    = Tgl          tangent of the track momentum dip angle
  // fP4    = C            track curvature

  // Magnetic field
  TVector3 field(0.,0.,1/AliKalmanTrack::GetConvConst());
  // radius [cm] of track projection in (x,y) 
  Double_t rho =
    mom.Pt()*AliKalmanTrack::GetConvConst();

  Double_t charge = 
    TDatabasePDG::Instance()->GetParticle(part->GetPdgCode())->Charge();
  charge /=3.;

  if (TMath::Abs(charge) < 0.9) charge=1.e-9;
 
  TVector3 vrho = mom.Cross(field);
  vrho *= charge;
  vrho.SetMag(rho);

  vrho += ver;
  Double_t x0 = vrho.X();

  xx[0] = ver.Y();
  xx[1] = ver.Z();
  xx[3] = mom.Pz()/mom.Pt();
  xx[4] = -charge/rho;
  xx[2] = xx[4]*x0;
  //  if (TMath::Abs(charge) < 0.9) xx[4] = 0;

  AliTPCtrack * tr = new AliTPCtrack(index, xx, cc, xref, alpha);
  tr->SetLabel(index);
  return tr;

}

void mcminiesd() {

  printf("====================\n");
  printf("Main program\n");
  TStopwatch timer;
  timer.Start();

  // Input "data" file, tree, and branch
  TFile * fIn = TFile::Open("AliESDs.root");
  TTree * tIn = (TTree*)fIn->Get("esdTree");
  TBranch * bIn = tIn->GetBranch("ESD");  

  AliESD * esdIn = 0; // The input ESD object is put here
  bIn->SetAddress(&esdIn);

  // Output file, trees, and branches.
  // Both "data" and MC are stored in the same file
  TFile * fOut   = TFile::Open("miniesd.root","recreate");
  TTree * tOut   = new TTree("esdTree","esdTree");
  TTree * tMC = new TTree("mcStackTree","mcStackTree");

  AliESD * esdOut = 0; // The newly created ESD object
  TBranch * bOut = tOut->Branch("ESD","AliESD",&esdOut);
  bOut->SetCompressionLevel(9);

  // TParticles are stored in TClonesArray
  // The corresponding branch is created using "dummy" TClonesArray
  TClonesArray * ministack = new TClonesArray("TParticle");
  TBranch * bMC = tMC->Branch("Stack",&ministack);
  ministack->Delete();
  delete ministack;
  bMC->SetCompressionLevel(9);

  // Event header
  AliHeader * headerIn = 0; //
  TBranch * bHeader = tMC->Branch("Header","AliHeader",&headerIn);
  bHeader->SetCompressionLevel(9);

 // Run loader
  AliRunLoader* runLoader = AliRunLoader::Open("galice.root");

  // gAlice
  runLoader->LoadgAlice();
  gAlice = runLoader->GetAliRun();

  // Magnetic field
  AliTracker::SetFieldMap(gAlice->Field(),1); // 1 means uniform magnetic field

  // Now load kinematics and event header
  runLoader->LoadKinematics();
  runLoader->LoadHeader();

  // Loop on events: check that MC and data contain the same number of events
  Int_t nevESD = bIn->GetEntries();
  Int_t nevMC = (Int_t)runLoader->GetNumberOfEvents();

  if (nevESD != nevMC)
    cout << "Inconsistent number of ESD("<<nevESD<<") and MC("<<nevMC<<") events" << endl;

  // Loop on events
  Int_t nev=TMath::Min(nevMC,nevESD);

  cout << nev << " events to process" << endl;

  for (Int_t iev=0; iev<nev; iev++) {
    cout << "Event " << iev << endl;
    // Get "data" event
    bIn->GetEntry(iev);

    // "Data" selection
    selectMiniESD(esdIn,esdOut);

    // Get MC event
    runLoader->GetEvent(iev);

    // MC selection
    kinetree(runLoader,esdOut,ministack);
    bMC->SetAddress(&ministack); // needed because ministack has changed

    // Header
    headerIn = runLoader->GetHeader();

    // Fill the trees
    tOut->Fill();
    tMC->Fill();
    delete esdOut;
    esdOut = 0x0;
    //    delete esdIn;
    //    esdIn = 0x0;
    ministack->Delete();
    delete ministack;

  }

  fOut = tOut->GetCurrentFile();
  fOut->Write();
  fOut->Close();

  fIn->Close();

  timer.Stop();
  timer.Print();
  printf("====================\n");
}
 mcminiesd.C:1
 mcminiesd.C:2
 mcminiesd.C:3
 mcminiesd.C:4
 mcminiesd.C:5
 mcminiesd.C:6
 mcminiesd.C:7
 mcminiesd.C:8
 mcminiesd.C:9
 mcminiesd.C:10
 mcminiesd.C:11
 mcminiesd.C:12
 mcminiesd.C:13
 mcminiesd.C:14
 mcminiesd.C:15
 mcminiesd.C:16
 mcminiesd.C:17
 mcminiesd.C:18
 mcminiesd.C:19
 mcminiesd.C:20
 mcminiesd.C:21
 mcminiesd.C:22
 mcminiesd.C:23
 mcminiesd.C:24
 mcminiesd.C:25
 mcminiesd.C:26
 mcminiesd.C:27
 mcminiesd.C:28
 mcminiesd.C:29
 mcminiesd.C:30
 mcminiesd.C:31
 mcminiesd.C:32
 mcminiesd.C:33
 mcminiesd.C:34
 mcminiesd.C:35
 mcminiesd.C:36
 mcminiesd.C:37
 mcminiesd.C:38
 mcminiesd.C:39
 mcminiesd.C:40
 mcminiesd.C:41
 mcminiesd.C:42
 mcminiesd.C:43
 mcminiesd.C:44
 mcminiesd.C:45
 mcminiesd.C:46
 mcminiesd.C:47
 mcminiesd.C:48
 mcminiesd.C:49
 mcminiesd.C:50
 mcminiesd.C:51
 mcminiesd.C:52
 mcminiesd.C:53
 mcminiesd.C:54
 mcminiesd.C:55
 mcminiesd.C:56
 mcminiesd.C:57
 mcminiesd.C:58
 mcminiesd.C:59
 mcminiesd.C:60
 mcminiesd.C:61
 mcminiesd.C:62
 mcminiesd.C:63
 mcminiesd.C:64
 mcminiesd.C:65
 mcminiesd.C:66
 mcminiesd.C:67
 mcminiesd.C:68
 mcminiesd.C:69
 mcminiesd.C:70
 mcminiesd.C:71
 mcminiesd.C:72
 mcminiesd.C:73
 mcminiesd.C:74
 mcminiesd.C:75
 mcminiesd.C:76
 mcminiesd.C:77
 mcminiesd.C:78
 mcminiesd.C:79
 mcminiesd.C:80
 mcminiesd.C:81
 mcminiesd.C:82
 mcminiesd.C:83
 mcminiesd.C:84
 mcminiesd.C:85
 mcminiesd.C:86
 mcminiesd.C:87
 mcminiesd.C:88
 mcminiesd.C:89
 mcminiesd.C:90
 mcminiesd.C:91
 mcminiesd.C:92
 mcminiesd.C:93
 mcminiesd.C:94
 mcminiesd.C:95
 mcminiesd.C:96
 mcminiesd.C:97
 mcminiesd.C:98
 mcminiesd.C:99
 mcminiesd.C:100
 mcminiesd.C:101
 mcminiesd.C:102
 mcminiesd.C:103
 mcminiesd.C:104
 mcminiesd.C:105
 mcminiesd.C:106
 mcminiesd.C:107
 mcminiesd.C:108
 mcminiesd.C:109
 mcminiesd.C:110
 mcminiesd.C:111
 mcminiesd.C:112
 mcminiesd.C:113
 mcminiesd.C:114
 mcminiesd.C:115
 mcminiesd.C:116
 mcminiesd.C:117
 mcminiesd.C:118
 mcminiesd.C:119
 mcminiesd.C:120
 mcminiesd.C:121
 mcminiesd.C:122
 mcminiesd.C:123
 mcminiesd.C:124
 mcminiesd.C:125
 mcminiesd.C:126
 mcminiesd.C:127
 mcminiesd.C:128
 mcminiesd.C:129
 mcminiesd.C:130
 mcminiesd.C:131
 mcminiesd.C:132
 mcminiesd.C:133
 mcminiesd.C:134
 mcminiesd.C:135
 mcminiesd.C:136
 mcminiesd.C:137
 mcminiesd.C:138
 mcminiesd.C:139
 mcminiesd.C:140
 mcminiesd.C:141
 mcminiesd.C:142
 mcminiesd.C:143
 mcminiesd.C:144
 mcminiesd.C:145
 mcminiesd.C:146
 mcminiesd.C:147
 mcminiesd.C:148
 mcminiesd.C:149
 mcminiesd.C:150
 mcminiesd.C:151
 mcminiesd.C:152
 mcminiesd.C:153
 mcminiesd.C:154
 mcminiesd.C:155
 mcminiesd.C:156
 mcminiesd.C:157
 mcminiesd.C:158
 mcminiesd.C:159
 mcminiesd.C:160
 mcminiesd.C:161
 mcminiesd.C:162
 mcminiesd.C:163
 mcminiesd.C:164
 mcminiesd.C:165
 mcminiesd.C:166
 mcminiesd.C:167
 mcminiesd.C:168
 mcminiesd.C:169
 mcminiesd.C:170
 mcminiesd.C:171
 mcminiesd.C:172
 mcminiesd.C:173
 mcminiesd.C:174
 mcminiesd.C:175
 mcminiesd.C:176
 mcminiesd.C:177
 mcminiesd.C:178
 mcminiesd.C:179
 mcminiesd.C:180
 mcminiesd.C:181
 mcminiesd.C:182
 mcminiesd.C:183
 mcminiesd.C:184
 mcminiesd.C:185
 mcminiesd.C:186
 mcminiesd.C:187
 mcminiesd.C:188
 mcminiesd.C:189
 mcminiesd.C:190
 mcminiesd.C:191
 mcminiesd.C:192
 mcminiesd.C:193
 mcminiesd.C:194
 mcminiesd.C:195
 mcminiesd.C:196
 mcminiesd.C:197
 mcminiesd.C:198
 mcminiesd.C:199
 mcminiesd.C:200
 mcminiesd.C:201
 mcminiesd.C:202
 mcminiesd.C:203
 mcminiesd.C:204
 mcminiesd.C:205
 mcminiesd.C:206
 mcminiesd.C:207
 mcminiesd.C:208
 mcminiesd.C:209
 mcminiesd.C:210
 mcminiesd.C:211
 mcminiesd.C:212
 mcminiesd.C:213
 mcminiesd.C:214
 mcminiesd.C:215
 mcminiesd.C:216
 mcminiesd.C:217
 mcminiesd.C:218
 mcminiesd.C:219
 mcminiesd.C:220
 mcminiesd.C:221
 mcminiesd.C:222
 mcminiesd.C:223
 mcminiesd.C:224
 mcminiesd.C:225
 mcminiesd.C:226
 mcminiesd.C:227
 mcminiesd.C:228
 mcminiesd.C:229
 mcminiesd.C:230
 mcminiesd.C:231
 mcminiesd.C:232
 mcminiesd.C:233
 mcminiesd.C:234
 mcminiesd.C:235
 mcminiesd.C:236
 mcminiesd.C:237
 mcminiesd.C:238
 mcminiesd.C:239
 mcminiesd.C:240
 mcminiesd.C:241
 mcminiesd.C:242
 mcminiesd.C:243
 mcminiesd.C:244
 mcminiesd.C:245
 mcminiesd.C:246
 mcminiesd.C:247
 mcminiesd.C:248
 mcminiesd.C:249
 mcminiesd.C:250
 mcminiesd.C:251
 mcminiesd.C:252
 mcminiesd.C:253
 mcminiesd.C:254
 mcminiesd.C:255
 mcminiesd.C:256
 mcminiesd.C:257
 mcminiesd.C:258
 mcminiesd.C:259
 mcminiesd.C:260
 mcminiesd.C:261
 mcminiesd.C:262
 mcminiesd.C:263
 mcminiesd.C:264
 mcminiesd.C:265
 mcminiesd.C:266
 mcminiesd.C:267
 mcminiesd.C:268
 mcminiesd.C:269
 mcminiesd.C:270
 mcminiesd.C:271
 mcminiesd.C:272
 mcminiesd.C:273
 mcminiesd.C:274
 mcminiesd.C:275
 mcminiesd.C:276
 mcminiesd.C:277
 mcminiesd.C:278
 mcminiesd.C:279
 mcminiesd.C:280
 mcminiesd.C:281
 mcminiesd.C:282
 mcminiesd.C:283
 mcminiesd.C:284
 mcminiesd.C:285
 mcminiesd.C:286
 mcminiesd.C:287
 mcminiesd.C:288
 mcminiesd.C:289
 mcminiesd.C:290
 mcminiesd.C:291
 mcminiesd.C:292
 mcminiesd.C:293
 mcminiesd.C:294
 mcminiesd.C:295
 mcminiesd.C:296
 mcminiesd.C:297
 mcminiesd.C:298
 mcminiesd.C:299
 mcminiesd.C:300
 mcminiesd.C:301
 mcminiesd.C:302
 mcminiesd.C:303
 mcminiesd.C:304
 mcminiesd.C:305
 mcminiesd.C:306
 mcminiesd.C:307
 mcminiesd.C:308
 mcminiesd.C:309
 mcminiesd.C:310
 mcminiesd.C:311
 mcminiesd.C:312
 mcminiesd.C:313
 mcminiesd.C:314
 mcminiesd.C:315
 mcminiesd.C:316
 mcminiesd.C:317
 mcminiesd.C:318
 mcminiesd.C:319
 mcminiesd.C:320
 mcminiesd.C:321
 mcminiesd.C:322
 mcminiesd.C:323
 mcminiesd.C:324
 mcminiesd.C:325
 mcminiesd.C:326
 mcminiesd.C:327
 mcminiesd.C:328
 mcminiesd.C:329
 mcminiesd.C:330
 mcminiesd.C:331
 mcminiesd.C:332
 mcminiesd.C:333
 mcminiesd.C:334
 mcminiesd.C:335
 mcminiesd.C:336
 mcminiesd.C:337
 mcminiesd.C:338
 mcminiesd.C:339
 mcminiesd.C:340
 mcminiesd.C:341
 mcminiesd.C:342
 mcminiesd.C:343
 mcminiesd.C:344
 mcminiesd.C:345
 mcminiesd.C:346
 mcminiesd.C:347
 mcminiesd.C:348
 mcminiesd.C:349
 mcminiesd.C:350
 mcminiesd.C:351
 mcminiesd.C:352
 mcminiesd.C:353
 mcminiesd.C:354
 mcminiesd.C:355
 mcminiesd.C:356
 mcminiesd.C:357
 mcminiesd.C:358
 mcminiesd.C:359
 mcminiesd.C:360
 mcminiesd.C:361
 mcminiesd.C:362
 mcminiesd.C:363
 mcminiesd.C:364
 mcminiesd.C:365
 mcminiesd.C:366
 mcminiesd.C:367
 mcminiesd.C:368
 mcminiesd.C:369
 mcminiesd.C:370
 mcminiesd.C:371
 mcminiesd.C:372
 mcminiesd.C:373
 mcminiesd.C:374
 mcminiesd.C:375
 mcminiesd.C:376
 mcminiesd.C:377
 mcminiesd.C:378
 mcminiesd.C:379
 mcminiesd.C:380
 mcminiesd.C:381
 mcminiesd.C:382
 mcminiesd.C:383
 mcminiesd.C:384
 mcminiesd.C:385
 mcminiesd.C:386
 mcminiesd.C:387
 mcminiesd.C:388
 mcminiesd.C:389
 mcminiesd.C:390
 mcminiesd.C:391
 mcminiesd.C:392
 mcminiesd.C:393
 mcminiesd.C:394
 mcminiesd.C:395
 mcminiesd.C:396
 mcminiesd.C:397
 mcminiesd.C:398
 mcminiesd.C:399
 mcminiesd.C:400
 mcminiesd.C:401
 mcminiesd.C:402
 mcminiesd.C:403
 mcminiesd.C:404
 mcminiesd.C:405
 mcminiesd.C:406
 mcminiesd.C:407
 mcminiesd.C:408
 mcminiesd.C:409
 mcminiesd.C:410
 mcminiesd.C:411
 mcminiesd.C:412
 mcminiesd.C:413
 mcminiesd.C:414
 mcminiesd.C:415
 mcminiesd.C:416
 mcminiesd.C:417
 mcminiesd.C:418
 mcminiesd.C:419
 mcminiesd.C:420
 mcminiesd.C:421
 mcminiesd.C:422
 mcminiesd.C:423
 mcminiesd.C:424
 mcminiesd.C:425
 mcminiesd.C:426
 mcminiesd.C:427
 mcminiesd.C:428
 mcminiesd.C:429
 mcminiesd.C:430
 mcminiesd.C:431
 mcminiesd.C:432
 mcminiesd.C:433
 mcminiesd.C:434
 mcminiesd.C:435
 mcminiesd.C:436
 mcminiesd.C:437
 mcminiesd.C:438
 mcminiesd.C:439
 mcminiesd.C:440
 mcminiesd.C:441
 mcminiesd.C:442
 mcminiesd.C:443
 mcminiesd.C:444
 mcminiesd.C:445
 mcminiesd.C:446
 mcminiesd.C:447
 mcminiesd.C:448
 mcminiesd.C:449
 mcminiesd.C:450
 mcminiesd.C:451
 mcminiesd.C:452
 mcminiesd.C:453
 mcminiesd.C:454
 mcminiesd.C:455
 mcminiesd.C:456
 mcminiesd.C:457
 mcminiesd.C:458
 mcminiesd.C:459
 mcminiesd.C:460
 mcminiesd.C:461
 mcminiesd.C:462
 mcminiesd.C:463
 mcminiesd.C:464
 mcminiesd.C:465
 mcminiesd.C:466
 mcminiesd.C:467
 mcminiesd.C:468
 mcminiesd.C:469
 mcminiesd.C:470
 mcminiesd.C:471
 mcminiesd.C:472
 mcminiesd.C:473
 mcminiesd.C:474
 mcminiesd.C:475
 mcminiesd.C:476
 mcminiesd.C:477
 mcminiesd.C:478
 mcminiesd.C:479
 mcminiesd.C:480
 mcminiesd.C:481
 mcminiesd.C:482
 mcminiesd.C:483
 mcminiesd.C:484
 mcminiesd.C:485
 mcminiesd.C:486
 mcminiesd.C:487
 mcminiesd.C:488
 mcminiesd.C:489
 mcminiesd.C:490
 mcminiesd.C:491
 mcminiesd.C:492
 mcminiesd.C:493
 mcminiesd.C:494
 mcminiesd.C:495
 mcminiesd.C:496
 mcminiesd.C:497
 mcminiesd.C:498
 mcminiesd.C:499
 mcminiesd.C:500
 mcminiesd.C:501
 mcminiesd.C:502
 mcminiesd.C:503
 mcminiesd.C:504
 mcminiesd.C:505
 mcminiesd.C:506
 mcminiesd.C:507
 mcminiesd.C:508
 mcminiesd.C:509
 mcminiesd.C:510
 mcminiesd.C:511
 mcminiesd.C:512
 mcminiesd.C:513
 mcminiesd.C:514
 mcminiesd.C:515
 mcminiesd.C:516
 mcminiesd.C:517
 mcminiesd.C:518
 mcminiesd.C:519
 mcminiesd.C:520
 mcminiesd.C:521
 mcminiesd.C:522
 mcminiesd.C:523
 mcminiesd.C:524
 mcminiesd.C:525
 mcminiesd.C:526
 mcminiesd.C:527
 mcminiesd.C:528
 mcminiesd.C:529
 mcminiesd.C:530
 mcminiesd.C:531
 mcminiesd.C:532
 mcminiesd.C:533
 mcminiesd.C:534
 mcminiesd.C:535
 mcminiesd.C:536
 mcminiesd.C:537
 mcminiesd.C:538
 mcminiesd.C:539
 mcminiesd.C:540
 mcminiesd.C:541
 mcminiesd.C:542
 mcminiesd.C:543
 mcminiesd.C:544
 mcminiesd.C:545
 mcminiesd.C:546
 mcminiesd.C:547
 mcminiesd.C:548
 mcminiesd.C:549
 mcminiesd.C:550
 mcminiesd.C:551
 mcminiesd.C:552
 mcminiesd.C:553
 mcminiesd.C:554
 mcminiesd.C:555
 mcminiesd.C:556
 mcminiesd.C:557
 mcminiesd.C:558
 mcminiesd.C:559
 mcminiesd.C:560
 mcminiesd.C:561
 mcminiesd.C:562
 mcminiesd.C:563
 mcminiesd.C:564
 mcminiesd.C:565
 mcminiesd.C:566
 mcminiesd.C:567
 mcminiesd.C:568
 mcminiesd.C:569
 mcminiesd.C:570
 mcminiesd.C:571
 mcminiesd.C:572
 mcminiesd.C:573
 mcminiesd.C:574
 mcminiesd.C:575
 mcminiesd.C:576
 mcminiesd.C:577
 mcminiesd.C:578
 mcminiesd.C:579
 mcminiesd.C:580
 mcminiesd.C:581
 mcminiesd.C:582
 mcminiesd.C:583
 mcminiesd.C:584
 mcminiesd.C:585
 mcminiesd.C:586
 mcminiesd.C:587
 mcminiesd.C:588
 mcminiesd.C:589
 mcminiesd.C:590
 mcminiesd.C:591
 mcminiesd.C:592
 mcminiesd.C:593
 mcminiesd.C:594
 mcminiesd.C:595
 mcminiesd.C:596
 mcminiesd.C:597
 mcminiesd.C:598