Fork me on GitHub

Changeset 4006893 in git


Ignore:
Timestamp:
Jun 1, 2021, 8:10:27 PM (4 years ago)
Author:
Pavel Demin <pavel.demin@…>
Branches:
master
Children:
fc6bce2
Parents:
06630ef
Message:

fix order of decay products in DelphesHepMC3Reader

Location:
classes
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • classes/DelphesHepMC3Reader.cc

    r06630ef r4006893  
    5656DelphesHepMC3Reader::DelphesHepMC3Reader() :
    5757  fInputFile(0), fBuffer(0), fPDG(0),
    58   fVertexCounter(-1), fParticleCounter(-1)
     58  fVertexCounter(-2), fParticleCounter(-1)
    5959{
    6060  fBuffer = new char[kBufferSize];
     
    8181void DelphesHepMC3Reader::Clear()
    8282{
    83   fWeight.clear();
     83  fWeights.clear();
    8484  fMomentumCoefficient = 1.0;
    8585  fPositionCoefficient = 1.0;
    86   fVertexCounter = -1;
     86  fVertexCounter = -2;
    8787  fParticleCounter = -1;
    88   fVertexMap.clear();
     88  fVertices.clear();
     89  fParticles.clear();
     90  fInVertexMap.clear();
     91  fOutVertexMap.clear();
     92  fMotherMap.clear();
    8993  fDaughterMap.clear();
    9094}
     
    9498bool DelphesHepMC3Reader::EventReady()
    9599{
    96   return (fParticleCounter == 0);
     100  return (fVertexCounter == -1) && (fParticleCounter == 0);
    97101}
    98102
     
    119123    Clear();
    120124
    121     fX = 0.0;
    122     fY = 0.0;
    123     fZ = 0.0;
    124     fT = 0.0;
    125 
    126125    rc = bufferStream.ReadInt(fEventNumber)
    127126      && bufferStream.ReadInt(fVertexCounter)
     
    168167    while(bufferStream.ReadDbl(weight))
    169168    {
    170       fWeight.push_back(weight);
     169      fWeights.push_back(weight);
    171170    }
    172171  }
     
    257256  else if(key == 'V')
    258257  {
     258    fParticles.clear();
     259
    259260    fX = 0.0;
    260261    fY = 0.0;
     
    262263    fT = 0.0;
    263264
    264     fM1 = 0;
    265     fM2 = 0;
    266 
    267265    rc = bufferStream.ReadInt(fVertexCode)
    268266      && bufferStream.ReadInt(fVertexStatus);
     
    286284    while(bufferStream.ReadInt(code))
    287285    {
    288       if(code < fM1 || fM1 == 0) fM1 = code;
    289       if(code > fM2) fM2 = code;
    290       fVertexMap[code] = fVertexCode;
    291     }
    292 
    293     if(fM1 == fM2) fM2 = 0;
     286      fParticles.push_back(code);
     287      bufferStream.FindChr(',');
     288    }
    294289
    295290    if(bufferStream.FindChr('@'))
     
    307302      }
    308303    }
     304
     305    AnalyzeVertex(factory, fVertexCode);
    309306  }
    310307  else if(key == 'P' && fParticleCounter > 0)
     
    329326    }
    330327
    331     itDaughterMap = fDaughterMap.find(fOutVertexCode);
    332     if(itDaughterMap == fDaughterMap.end())
    333     {
    334       fDaughterMap[fOutVertexCode] = make_pair(fParticleCode, fParticleCode);
    335     }
    336     else
    337     {
    338       itDaughterMap->second.second = fParticleCode;
    339     }
    340 
    341328    AnalyzeParticle(factory, allParticleOutputArray,
    342329      stableParticleOutputArray, partonOutputArray);
     
    363350  element->ProcessID = fProcessID;
    364351  element->MPI = fMPI;
    365   element->Weight = fWeight.size() > 0 ? fWeight[0] : 1.0;
     352  element->Weight = fWeights.size() > 0 ? fWeights[0] : 1.0;
    366353  element->CrossSection = fCrossSection;
    367354  element->CrossSectionError = fCrossSectionError;
     
    389376  vector<double>::const_iterator itWeight;
    390377
    391   for(itWeight = fWeight.begin(); itWeight != fWeight.end(); ++itWeight)
     378  for(itWeight = fWeights.begin(); itWeight != fWeights.end(); ++itWeight)
    392379  {
    393380    element = static_cast<Weight *>(branch->NewEntry());
    394381
    395382    element->Weight = *itWeight;
     383  }
     384}
     385
     386//---------------------------------------------------------------------------
     387
     388void DelphesHepMC3Reader::AnalyzeVertex(DelphesFactory *factory, int code, Candidate *candidate)
     389{
     390  int index;
     391  TLorentzVector *position;
     392  TObjArray *array;
     393  vector<int>::iterator itParticle;
     394  map<int, int>::iterator itVertexMap;
     395
     396  itVertexMap = fOutVertexMap.find(code);
     397  if(itVertexMap == fOutVertexMap.end())
     398  {
     399    --fVertexCounter;
     400
     401    index = fVertices.size();
     402    fOutVertexMap[code] = index;
     403    if(candidate && code > 0) fInVertexMap[code] = index;
     404
     405    position = factory->New<TLorentzVector>();
     406    array = factory->NewArray();
     407    position->SetXYZT(0.0, 0.0, 0.0, 0.0);
     408    array->Clear();
     409    fVertices.push_back(make_pair(position, array));
     410  }
     411  else
     412  {
     413    index = itVertexMap->second;
     414    position = fVertices[index].first;
     415    array = fVertices[index].second;
     416  }
     417
     418  if(candidate)
     419  {
     420    array->Add(candidate);
     421  }
     422  else
     423  {
     424    position->SetXYZT(fX, fY, fZ, fT);
     425    for(itParticle = fParticles.begin(); itParticle != fParticles.end(); ++itParticle)
     426    {
     427      fInVertexMap[*itParticle] = index;
     428    }
    396429  }
    397430}
     
    425458  }
    426459
    427   candidate->Position.SetXYZT(fX, fY, fZ, fT);
    428   if(fPositionCoefficient != 1.0)
    429   {
    430     candidate->Position *= fPositionCoefficient;
    431   }
    432 
    433   candidate->D1 = -1;
    434   candidate->D2 = -1;
    435 
    436   if(fOutVertexCode < 0)
    437   {
    438     candidate->M1 = fM1 - 1;
    439     candidate->M2 = fM2 - 1;
    440   }
    441   else
    442   {
    443     candidate->M1 = fOutVertexCode - 1;
    444     candidate->M2 = -1;
    445   }
    446 
    447   allParticleOutputArray->Add(candidate);
     460  candidate->D1 = fParticleCode;
     461
     462  AnalyzeVertex(factory, fOutVertexCode, candidate);
    448463
    449464  if(!pdgParticle) return;
     
    463478void DelphesHepMC3Reader::FinalizeParticles(TObjArray *allParticleOutputArray)
    464479{
     480  TLorentzVector *position;
     481  TObjArray *array;
    465482  Candidate *candidate;
    466483  map<int, int >::iterator itVertexMap;
     484  map<int, pair<int, int> >::iterator itMotherMap;
    467485  map<int, pair<int, int> >::iterator itDaughterMap;
    468   int i, index;
     486  int i, j, code, counter;
     487
     488  counter = 0;
     489  for(i = 0; i < fVertices.size(); ++i)
     490  {
     491    position = fVertices[i].first;
     492    array = fVertices[i].second;
     493
     494    for(j = 0; j < array->GetEntriesFast(); ++j)
     495    {
     496      candidate = static_cast<Candidate *>(array->At(j));
     497
     498      candidate->Position = *position;
     499      if(fPositionCoefficient != 1.0)
     500      {
     501        candidate->Position *= fPositionCoefficient;
     502      }
     503
     504      candidate->M1 = i;
     505
     506      itDaughterMap = fDaughterMap.find(i);
     507      if(itDaughterMap == fDaughterMap.end())
     508      {
     509        fDaughterMap[i] = make_pair(counter, counter);
     510      }
     511      else
     512      {
     513        itDaughterMap->second.second = counter;
     514      }
     515
     516      code = candidate->D1;
     517
     518      itVertexMap = fInVertexMap.find(code);
     519      if(itVertexMap == fInVertexMap.end())
     520      {
     521        candidate->D1 = -1;
     522      }
     523      else
     524      {
     525        code = itVertexMap->second;
     526
     527        candidate->D1 = code;
     528
     529        itMotherMap = fMotherMap.find(code);
     530        if(itMotherMap == fMotherMap.end())
     531        {
     532          fMotherMap[code] = make_pair(counter, -1);
     533        }
     534        else
     535        {
     536          itMotherMap->second.second = counter;
     537        }
     538      }
     539
     540      allParticleOutputArray->Add(candidate);
     541
     542      ++counter;
     543    }
     544  }
    469545
    470546  for(i = 0; i < allParticleOutputArray->GetEntriesFast(); ++i)
     
    472548    candidate = static_cast<Candidate *>(allParticleOutputArray->At(i));
    473549
    474     index = i + 1;
    475 
    476     itVertexMap = fVertexMap.find(index);
    477     if(itVertexMap != fVertexMap.end())
    478     {
    479       index = itVertexMap->second;
    480     }
    481 
    482     itDaughterMap = fDaughterMap.find(index);
    483     if(itDaughterMap == fDaughterMap.end())
     550    itMotherMap = fMotherMap.find(candidate->M1);
     551    if(itMotherMap == fMotherMap.end())
     552    {
     553      candidate->M1 = -1;
     554      candidate->M2 = -1;
     555    }
     556    else
     557    {
     558      candidate->M1 = itMotherMap->second.first;
     559      candidate->M2 = itMotherMap->second.second;
     560    }
     561
     562    if(candidate->D1 < 0)
    484563    {
    485564      candidate->D1 = -1;
     
    488567    else
    489568    {
    490       candidate->D1 = itDaughterMap->second.first - 1;
    491       candidate->D2 = itDaughterMap->second.second - 1;
    492     }
    493   }
    494 }
    495 
    496 //---------------------------------------------------------------------------
     569      itDaughterMap = fDaughterMap.find(candidate->D1);
     570      if(itDaughterMap == fDaughterMap.end())
     571      {
     572        candidate->D1 = -1;
     573        candidate->D2 = -1;
     574      }
     575      else
     576      {
     577        candidate->D1 = itDaughterMap->second.first;
     578        candidate->D2 = itDaughterMap->second.second;
     579      }
     580    }
     581  }
     582}
     583
     584//---------------------------------------------------------------------------
  • classes/DelphesHepMC3Reader.h

    r06630ef r4006893  
    3636class TStopwatch;
    3737class TDatabasePDG;
     38class TLorentzVector;
    3839class ExRootTreeBranch;
    3940class DelphesFactory;
     41class Candidate;
    4042
    4143class DelphesHepMC3Reader
     
    6163
    6264private:
     65  void AnalyzeVertex(DelphesFactory *factory, int code, Candidate *candidate = 0);
     66
    6367  void AnalyzeParticle(DelphesFactory *factory,
    6468    TObjArray *allParticleOutputArray,
     
    7983  double fMomentumCoefficient, fPositionCoefficient;
    8084
    81   std::vector<double> fWeight;
     85  std::vector<double> fWeights;
    8286
    8387  double fCrossSection, fCrossSectionError;
     
    9296  double fPx, fPy, fPz, fE, fMass;
    9397
    94   int fM1, fM2;
     98  std::vector<std::pair<TLorentzVector *, TObjArray *> > fVertices;
     99  std::vector<int> fParticles;
    95100
    96   std::map<int, int> fVertexMap;
     101  std::map<int, int> fInVertexMap;
     102  std::map<int, int> fOutVertexMap;
     103
     104  std::map<int, std::pair<int, int> > fMotherMap;
    97105  std::map<int, std::pair<int, int> > fDaughterMap;
    98106};
Note: See TracChangeset for help on using the changeset viewer.