Fork me on GitHub

Ignore:
File:
1 edited

Legend:

Unmodified
Added
Removed
  • modules/Calorimeter.cc

    rf8299bc r298734e  
    5858  fItParticleInputArray(0), fItTrackInputArray(0)
    5959{
    60   Int_t i;
    61 
     60 
    6261  fECalResolutionFormula = new DelphesFormula;
    6362  fHCalResolutionFormula = new DelphesFormula;
    6463
    65   for(i = 0; i < 2; ++i)
    66   {
    67     fECalTowerTrackArray[i] = new TObjArray;
    68     fItECalTowerTrackArray[i] = fECalTowerTrackArray[i]->MakeIterator();
    69 
    70     fHCalTowerTrackArray[i] = new TObjArray;
    71     fItHCalTowerTrackArray[i] = fHCalTowerTrackArray[i]->MakeIterator();
    72   }
     64  fECalTowerTrackArray = new TObjArray;
     65  fItECalTowerTrackArray = fECalTowerTrackArray->MakeIterator();
     66
     67  fHCalTowerTrackArray = new TObjArray;
     68  fItHCalTowerTrackArray = fHCalTowerTrackArray->MakeIterator();
     69 
    7370}
    7471
     
    7774Calorimeter::~Calorimeter()
    7875{
    79   Int_t i;
    80 
     76 
    8177  if(fECalResolutionFormula) delete fECalResolutionFormula;
    8278  if(fHCalResolutionFormula) delete fHCalResolutionFormula;
    8379
    84   for(i = 0; i < 2; ++i)
    85   {
    86     if(fECalTowerTrackArray[i]) delete fECalTowerTrackArray[i];
    87     if(fItECalTowerTrackArray[i]) delete fItECalTowerTrackArray[i];
    88 
    89     if(fHCalTowerTrackArray[i]) delete fHCalTowerTrackArray[i];
    90     if(fItHCalTowerTrackArray[i]) delete fItHCalTowerTrackArray[i];
    91   }
     80  if(fECalTowerTrackArray) delete fECalTowerTrackArray;
     81  if(fItECalTowerTrackArray) delete fItECalTowerTrackArray;
     82
     83  if(fHCalTowerTrackArray) delete fHCalTowerTrackArray;
     84  if(fItHCalTowerTrackArray) delete fItHCalTowerTrackArray;
     85 
    9286}
    9387
     
    219213  Double_t ecalEnergy, hcalEnergy;
    220214  Double_t ecalSigma, hcalSigma;
     215  Double_t energyGuess;
    221216  Int_t pdgCode;
    222217
     
    368363      fHCalTowerEnergy = 0.0;
    369364
    370       fECalTrackEnergy[0] = 0.0;
    371       fECalTrackEnergy[1] = 0.0;
    372 
    373       fHCalTrackEnergy[0] = 0.0;
    374       fHCalTrackEnergy[1] = 0.0;
    375 
     365      fECalTrackEnergy = 0.0;
     366      fHCalTrackEnergy = 0.0;
     367     
     368      fECalTrackSigma = 0.0;
     369      fHCalTrackSigma = 0.0;
     370     
    376371      fTowerTrackHits = 0;
    377372      fTowerPhotonHits = 0;
    378373
    379       fECalTowerTrackArray[0]->Clear();
    380       fECalTowerTrackArray[1]->Clear();
    381 
    382       fHCalTowerTrackArray[0]->Clear();
    383       fHCalTowerTrackArray[1]->Clear();
     374      fECalTowerTrackArray->Clear();
     375      fHCalTowerTrackArray->Clear();
     376   
    384377    }
    385378
     
    406399      if(fECalTrackFractions[number] > 1.0E-9 && fHCalTrackFractions[number] < 1.0E-9)
    407400      {
    408         ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
    409         if(ecalSigma/momentum.E() < track->TrackResolution)
    410         {
    411           fECalTrackEnergy[0] += ecalEnergy;
    412           fECalTowerTrackArray[0]->Add(track);
    413         }
    414         else
    415         {
    416           fECalTrackEnergy[1] += ecalEnergy;
    417           fECalTowerTrackArray[1]->Add(track);
    418         }
     401        fECalTrackEnergy += ecalEnergy;
     402        ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());       
     403        if(ecalSigma/momentum.E() < track->TrackResolution) energyGuess = ecalEnergy;       
     404        else energyGuess = momentum.E();
     405
     406        fECalTrackSigma += (track->TrackResolution)*energyGuess*(track->TrackResolution)*energyGuess;
     407        fECalTowerTrackArray->Add(track);
    419408      }
     409     
    420410      else if(fECalTrackFractions[number] < 1.0E-9 && fHCalTrackFractions[number] > 1.0E-9)
    421411      {
     412        fHCalTrackEnergy += hcalEnergy;
    422413        hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
    423         if(hcalSigma/momentum.E() < track->TrackResolution)
    424         {
    425           fHCalTrackEnergy[0] += hcalEnergy;
    426           fHCalTowerTrackArray[0]->Add(track);
    427         }
    428         else
    429         {
    430           fHCalTrackEnergy[1] += hcalEnergy;
    431           fHCalTowerTrackArray[1]->Add(track);
    432         }
     414        if(hcalSigma/momentum.E() < track->TrackResolution) energyGuess = hcalEnergy;
     415        else energyGuess = momentum.E();
     416
     417        fHCalTrackSigma += (track->TrackResolution)*energyGuess*(track->TrackResolution)*energyGuess;
     418        fHCalTowerTrackArray->Add(track);
    433419      }
     420     
    434421      else if(fECalTrackFractions[number] < 1.0E-9 && fHCalTrackFractions[number] < 1.0E-9)
    435422      {
     
    476463  Double_t energy, pt, eta, phi;
    477464  Double_t ecalEnergy, hcalEnergy;
     465  Double_t ecalNeutralEnergy, hcalNeutralEnergy;
     466 
    478467  Double_t ecalSigma, hcalSigma;
    479 
     468  Double_t ecalNeutralSigma, hcalNeutralSigma;
     469
     470  Double_t weightTrack, weightCalo, bestEnergyEstimate, rescaleFactor;
     471 
    480472  TLorentzVector momentum;
    481473  TFractionMap::iterator itFractionMap;
     
    484476
    485477  if(!fTower) return;
    486 
    487478
    488479  ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fECalTowerEnergy);
     
    554545    fTowerOutputArray->Add(fTower);
    555546  }
    556 
     547 
    557548  // fill energy flow candidates
    558 
    559   ecalEnergy -= fECalTrackEnergy[1];
    560   hcalEnergy -= fHCalTrackEnergy[1];
    561 
    562   fItECalTowerTrackArray[0]->Reset();
    563   while((track = static_cast<Candidate*>(fItECalTowerTrackArray[0]->Next())))
    564   {
    565     mother = track;
    566     track = static_cast<Candidate*>(track->Clone());
    567     track->AddCandidate(mother);
    568 
    569     track->Momentum *= ecalEnergy/fECalTrackEnergy[0];
    570 
    571     fEFlowTrackOutputArray->Add(track);
    572   }
    573 
    574   fItECalTowerTrackArray[1]->Reset();
    575   while((track = static_cast<Candidate*>(fItECalTowerTrackArray[1]->Next())))
    576   {
    577     mother = track;
    578     track = static_cast<Candidate*>(track->Clone());
    579     track->AddCandidate(mother);
    580 
    581     fEFlowTrackOutputArray->Add(track);
    582   }
    583 
    584   fItHCalTowerTrackArray[0]->Reset();
    585   while((track = static_cast<Candidate*>(fItHCalTowerTrackArray[0]->Next())))
    586   {
    587     mother = track;
    588     track = static_cast<Candidate*>(track->Clone());
    589     track->AddCandidate(mother);
    590 
    591     track->Momentum *= hcalEnergy/fHCalTrackEnergy[0];
    592 
    593     fEFlowTrackOutputArray->Add(track);
    594   }
    595 
    596   fItHCalTowerTrackArray[1]->Reset();
    597   while((track = static_cast<Candidate*>(fItHCalTowerTrackArray[1]->Next())))
    598   {
    599     mother = track;
    600     track = static_cast<Candidate*>(track->Clone());
    601     track->AddCandidate(mother);
    602 
    603     fEFlowTrackOutputArray->Add(track);
    604   }
    605 
    606   if(fECalTowerTrackArray[0]->GetEntriesFast() > 0) ecalEnergy = 0.0;
    607   if(fHCalTowerTrackArray[0]->GetEntriesFast() > 0) hcalEnergy = 0.0;
    608 
    609   ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, ecalEnergy);
    610   hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, hcalEnergy);
    611 
    612   if(ecalEnergy < fECalEnergyMin || ecalEnergy < fECalEnergySignificanceMin*ecalSigma) ecalEnergy = 0.0;
    613   if(hcalEnergy < fHCalEnergyMin || hcalEnergy < fHCalEnergySignificanceMin*hcalSigma) hcalEnergy = 0.0;
    614 
    615   energy = ecalEnergy + hcalEnergy;
    616 
    617   if(ecalEnergy > 0.0)
     549  fECalTrackSigma = TMath::Sqrt(fECalTrackSigma);
     550  fHCalTrackSigma = TMath::Sqrt(fHCalTrackSigma);
     551
     552  //compute neutral excesses
     553  ecalNeutralEnergy = max( (ecalEnergy - fECalTrackEnergy) , 0.0);
     554  hcalNeutralEnergy = max( (hcalEnergy - fHCalTrackEnergy) , 0.0);
     555 
     556  ecalNeutralSigma = ecalNeutralEnergy / TMath::Sqrt(fECalTrackSigma*fECalTrackSigma + ecalSigma*ecalSigma);
     557  hcalNeutralSigma = hcalNeutralEnergy / TMath::Sqrt(fHCalTrackSigma*fHCalTrackSigma + hcalSigma*hcalSigma);
     558 
     559   // if ecal neutral excess is significant, simply create neutral EflowPhoton tower and clone each track into eflowtrack
     560  if(ecalNeutralEnergy > fECalEnergyMin && ecalNeutralSigma > fECalEnergySignificanceMin)
    618561  {
    619562    // create new photon tower
    620563    tower = static_cast<Candidate*>(fTower->Clone());
    621 
    622     pt = ecalEnergy / TMath::CosH(eta);
    623 
    624     tower->Momentum.SetPtEtaPhiE(pt, eta, phi, ecalEnergy);
    625     tower->Eem = ecalEnergy;
     564    pt =  ecalNeutralEnergy / TMath::CosH(eta);
     565   
     566    tower->Momentum.SetPtEtaPhiE(pt, eta, phi, ecalNeutralEnergy);
     567    tower->Eem = ecalNeutralEnergy;
    626568    tower->Ehad = 0.0;
    627569    tower->PID = 22;
    628 
     570   
    629571    fEFlowPhotonOutputArray->Add(tower);
    630   }
    631   if(hcalEnergy > 0.0)
    632   {
    633     // create new neutral hadron tower
     572   
     573    //clone tracks
     574    fItECalTowerTrackArray->Reset();
     575    while((track = static_cast<Candidate*>(fItECalTowerTrackArray->Next())))
     576    {
     577      mother = track;
     578      track = static_cast<Candidate*>(track->Clone());
     579      track->AddCandidate(mother);
     580
     581      fEFlowTrackOutputArray->Add(track);
     582    }
     583 
     584  }
     585 
     586  // if neutral excess is not significant, rescale eflow tracks, such that the total charged equals the best measurement given by the calorimeter and tracking
     587  else if(fECalTrackEnergy > 0.0)
     588  {
     589    weightTrack = (fECalTrackSigma > 0.0) ? 1 / (fECalTrackSigma*fECalTrackSigma) : 0.0;
     590    weightCalo  = (ecalSigma > 0.0) ? 1 / (ecalSigma*ecalSigma) : 0.0;
     591 
     592    bestEnergyEstimate = (weightTrack*fECalTrackEnergy + weightCalo*ecalEnergy) / (weightTrack + weightCalo);
     593    rescaleFactor = bestEnergyEstimate/fECalTrackEnergy;
     594
     595    //rescale tracks
     596    fItECalTowerTrackArray->Reset();
     597    while((track = static_cast<Candidate*>(fItECalTowerTrackArray->Next())))
     598    { 
     599      mother = track;
     600      track = static_cast<Candidate*>(track->Clone());
     601      track->AddCandidate(mother);
     602
     603      track->Momentum *= rescaleFactor;
     604
     605      fEFlowTrackOutputArray->Add(track);
     606    }
     607  }
     608
     609
     610  // if hcal neutral excess is significant, simply create neutral EflowNeutralHadron tower and clone each track into eflowtrack
     611  if(hcalNeutralEnergy > fHCalEnergyMin && hcalNeutralSigma > fHCalEnergySignificanceMin)
     612  {
     613    // create new photon tower
    634614    tower = static_cast<Candidate*>(fTower->Clone());
    635 
    636     pt = hcalEnergy / TMath::CosH(eta);
    637 
    638     tower->Momentum.SetPtEtaPhiE(pt, eta, phi, hcalEnergy);
     615    pt =  hcalNeutralEnergy / TMath::CosH(eta);
     616   
     617    tower->Momentum.SetPtEtaPhiE(pt, eta, phi, hcalNeutralEnergy);
     618    tower->Ehad = hcalNeutralEnergy;
    639619    tower->Eem = 0.0;
    640     tower->Ehad = hcalEnergy;
    641 
     620   
    642621    fEFlowNeutralHadronOutputArray->Add(tower);
    643   }
     622   
     623    //clone tracks
     624    fItHCalTowerTrackArray->Reset();
     625    while((track = static_cast<Candidate*>(fItHCalTowerTrackArray->Next())))
     626    {
     627      mother = track;
     628      track = static_cast<Candidate*>(track->Clone());
     629      track->AddCandidate(mother);
     630
     631      fEFlowTrackOutputArray->Add(track);
     632    }
     633 
     634  }
     635 
     636  // if neutral excess is not significant, rescale eflow tracks, such that the total charged equals the best measurement given by the calorimeter and tracking
     637  else if(fHCalTrackEnergy > 0.0)
     638  {
     639    weightTrack = (fHCalTrackSigma > 0.0) ? 1 / (fHCalTrackSigma*fHCalTrackSigma) : 0.0;
     640    weightCalo  = (hcalSigma > 0.0) ? 1 / (hcalSigma*hcalSigma) : 0.0;
     641 
     642    bestEnergyEstimate = (weightTrack*fHCalTrackEnergy + weightCalo*hcalEnergy) / (weightTrack + weightCalo);
     643    rescaleFactor = bestEnergyEstimate / fHCalTrackEnergy;
     644
     645    //rescale tracks
     646    fItHCalTowerTrackArray->Reset();
     647    while((track = static_cast<Candidate*>(fItHCalTowerTrackArray->Next())))
     648    { 
     649      mother = track;
     650      track = static_cast<Candidate*>(track->Clone());
     651      track->AddCandidate(mother);
     652
     653      track->Momentum *= rescaleFactor;
     654
     655      fEFlowTrackOutputArray->Add(track);
     656    }
     657  }
     658 
     659 
    644660}
    645661
Note: See TracChangeset for help on using the changeset viewer.