Changeset d07e957 in git
- Timestamp:
- Nov 8, 2013, 5:09:58 PM (11 years ago)
- Branches:
- ImprovedOutputFile, Timing, dual_readout, llp, master
- Children:
- 87cadbd
- Parents:
- 984cd31
- Files:
-
- 7 edited
Legend:
- Unmodified
- Added
- Removed
-
classes/ClassesLinkDef.h
r984cd31 rd07e957 34 34 #pragma link C++ class HepMCEvent+; 35 35 #pragma link C++ class GenParticle+; 36 #pragma link C++ class Vertex+; 36 37 #pragma link C++ class MissingET+; 37 38 #pragma link C++ class ScalarHT+; -
classes/DelphesClasses.h
r984cd31 rd07e957 144 144 //--------------------------------------------------------------------------- 145 145 146 class Vertex: public TObject 147 { 148 public: 149 Float_t T; // vertex position (t component) 150 Float_t X; // vertex position (x component) 151 Float_t Y; // vertex position (y component) 152 Float_t Z; // vertex position (z component) 153 154 ClassDef(Vertex, 1) 155 }; 156 157 //--------------------------------------------------------------------------- 158 146 159 class MissingET: public TObject 147 160 { -
examples/delphes_card_CMS_PileUp.tcl
r984cd31 rd07e957 631 631 add Branch ScalarHT/energy ScalarHT ScalarHT 632 632 add Branch Rho/rho Rho Rho 633 } 634 633 add Branch PileUpMerger/Vertices Vertex Vertex 634 } 635 -
modules/PileUpMerger.cc
r984cd31 rd07e957 70 70 71 71 // create output arrays 72 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles")); 72 fParticleOutputArray = ExportArray(GetString("ParticleOutputArray", "stableParticles")); 73 fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices")); 73 74 } 74 75 … … 98 99 while((candidate = static_cast<Candidate*>(fItInputArray->Next()))) 99 100 { 100 f OutputArray->Add(candidate);101 fParticleOutputArray->Add(candidate); 101 102 } 102 103 … … 131 132 dphi = gRandom->Uniform(-TMath::Pi(), TMath::Pi()); 132 133 134 candidate = factory->NewCandidate(); 135 candidate->Position.SetXYZT(0.0, 0.0, dz, 0.0); 136 fVertexOutputArray->Add(candidate); 137 133 138 while(fReader->ReadParticle(pid, x, y, z, t, px, py, pz, e)) 134 139 { … … 151 156 candidate->Position.RotateZ(dphi); 152 157 153 f OutputArray->Add(candidate);158 fParticleOutputArray->Add(candidate); 154 159 } 155 160 } -
modules/PileUpMerger.h
r984cd31 rd07e957 43 43 const TObjArray *fInputArray; //! 44 44 45 TObjArray *fOutputArray; //! 45 TObjArray *fParticleOutputArray; //! 46 TObjArray *fVertexOutputArray; //! 46 47 47 48 ClassDef(PileUpMerger, 1) -
modules/TreeWriter.cc
r984cd31 rd07e957 213 213 //------------------------------------------------------------------------------ 214 214 215 void TreeWriter::ProcessVertices(ExRootTreeBranch *branch, TObjArray *array) 216 { 217 TIter iterator(array); 218 Candidate *candidate = 0; 219 Vertex *entry = 0; 220 221 // loop over all vertices 222 iterator.Reset(); 223 while((candidate = static_cast<Candidate*>(iterator.Next()))) 224 { 225 const TLorentzVector &position = candidate->Position; 226 227 entry = static_cast<Vertex*>(branch->NewEntry()); 228 229 entry->X = position.X(); 230 entry->Y = position.Y(); 231 entry->Z = position.Z(); 232 entry->T = position.T(); 233 } 234 } 235 236 //------------------------------------------------------------------------------ 237 215 238 void TreeWriter::ProcessTracks(ExRootTreeBranch *branch, TObjArray *array) 216 239 { … … 221 244 Double_t pt, signz, cosTheta, eta, rapidity; 222 245 223 // loop over all jets246 // loop over all tracks 224 247 iterator.Reset(); 225 248 while((candidate = static_cast<Candidate*>(iterator.Next()))) … … 280 303 Double_t pt, signPz, cosTheta, eta, rapidity; 281 304 282 // loop over all jets305 // loop over all towers 283 306 iterator.Reset(); 284 307 while((candidate = static_cast<Candidate*>(iterator.Next()))) -
modules/TreeWriter.h
r984cd31 rd07e957 41 41 42 42 void ProcessParticles(ExRootTreeBranch *branch, TObjArray *array); 43 void ProcessVertices(ExRootTreeBranch *branch, TObjArray *array); 43 44 void ProcessTracks(ExRootTreeBranch *branch, TObjArray *array); 44 45 void ProcessTowers(ExRootTreeBranch *branch, TObjArray *array);
Note:
See TracChangeset
for help on using the changeset viewer.