- Timestamp:
- Aug 16, 2013, 4:29:25 PM (11 years ago)
- Branches:
- ImprovedOutputFile, Timing, dual_readout, llp, master
- Children:
- 82e07cf
- Parents:
- 0a0a5ef
- Location:
- classes
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
classes/DelphesLHEFReader.cc
r0a0a5ef r0dc0eeb 41 41 DelphesLHEFReader::DelphesLHEFReader() : 42 42 fInputFile(0), fBuffer(0), fPDG(0), 43 fEvent Counter(-1), fParticleCounter(-1)43 fEventReady(kFALSE), fEventCounter(-1), fParticleCounter(-1) 44 44 { 45 45 fBuffer = new char[kBufferSize]; … … 66 66 void DelphesLHEFReader::Clear() 67 67 { 68 fEventReady = kFALSE; 68 69 fEventCounter = -1; 69 70 fParticleCounter = -1; 71 fRwgtList.clear(); 70 72 } 71 73 … … 74 76 bool DelphesLHEFReader::EventReady() 75 77 { 76 return (fEventCounter == 0) && (fParticleCounter == 0);78 return fEventReady; 77 79 } 78 80 … … 85 87 { 86 88 int rc; 89 char *pch; 90 double weight; 87 91 88 92 if(!fgets(fBuffer, kBufferSize, fInputFile)) return kFALSE; 89 90 DelphesStream bufferStream(fBuffer);91 93 92 94 if(strstr(fBuffer, "<event>")) … … 97 99 else if(fEventCounter > 0) 98 100 { 101 DelphesStream bufferStream(fBuffer); 102 99 103 rc = bufferStream.ReadInt(fParticleCounter) 100 104 && bufferStream.ReadInt(fProcessID) … … 114 118 else if(fParticleCounter > 0) 115 119 { 120 DelphesStream bufferStream(fBuffer); 121 116 122 rc = bufferStream.ReadInt(fPID) 117 123 && bufferStream.ReadInt(fStatus) … … 137 143 --fParticleCounter; 138 144 } 145 else if(strstr(fBuffer, "<wgt")) 146 { 147 pch = strstr(fBuffer, ">"); 148 if(!pch) 149 { 150 cerr << "** ERROR: " << "invalid weight format" << endl; 151 return kFALSE; 152 } 153 154 DelphesStream bufferStream(pch + 1); 155 rc = bufferStream.ReadDbl(weight); 156 157 if(!rc) 158 { 159 cerr << "** ERROR: " << "invalid weight format" << endl; 160 return kFALSE; 161 } 162 163 fRwgtList.push_back(weight); 164 } 165 else if(strstr(fBuffer, "</event>")) 166 { 167 fEventReady = kTRUE; 168 } 139 169 140 170 return kTRUE; … … 163 193 //--------------------------------------------------------------------------- 164 194 195 void DelphesLHEFReader::AnalyzeRwgt(ExRootTreeBranch *branch) 196 { 197 Weight *element; 198 vector<double>::const_iterator itRwgtList; 199 200 for(itRwgtList = fRwgtList.begin(); itRwgtList != fRwgtList.end(); ++itRwgtList) 201 { 202 element = static_cast<Weight *>(branch->NewEntry()); 203 204 element->Weight = *itRwgtList; 205 } 206 } 207 208 //--------------------------------------------------------------------------- 209 165 210 void DelphesLHEFReader::AnalyzeParticle(DelphesFactory *factory, 166 211 TObjArray *allParticleOutputArray, -
classes/DelphesLHEFReader.h
r0a0a5ef r0dc0eeb 16 16 17 17 #include <stdio.h> 18 19 #include <vector> 18 20 19 21 class TObjArray; … … 43 45 TStopwatch *readStopWatch, TStopwatch *procStopWatch); 44 46 47 void AnalyzeRwgt(ExRootTreeBranch *branch); 48 45 49 private: 46 50 … … 56 60 TDatabasePDG *fPDG; 57 61 62 bool fEventReady; 63 58 64 int fEventCounter; 59 65 … … 63 69 int fPID, fStatus, fM1, fM2, fC1, fC2; 64 70 double fPx, fPy, fPz, fE, fMass; 71 72 std::vector<double> fRwgtList; 65 73 }; 66 74
Note:
See TracChangeset
for help on using the changeset viewer.