Changeset 5496767 in git
- Timestamp:
- Jun 17, 2016, 6:52:57 PM (8 years ago)
- Branches:
- ImprovedOutputFile, Timing, dual_readout, llp, master
- Children:
- 62d3bc5
- Parents:
- 0f73bea
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
classes/DelphesClasses.h
r0f73bea r5496767 178 178 Float_t Y; // vertex position (y component) 179 179 Float_t Z; // vertex position (z component) 180 180 181 181 Double_t ErrorX; 182 182 Double_t ErrorY; … … 191 191 Double_t GenSumPT2; 192 192 193 TRefArray Constituents; // references to constituents 194 193 195 static CompBase *fgCompare; //! 194 196 const CompBase *GetCompare() const { return fgCompare; } … … 416 418 417 419 Float_t Eta; // track pseudorapidity 418 420 419 421 Float_t EtaOuter; // track pseudorapidity at the tracker edge 420 422 Float_t PhiOuter; // track azimuthal angle at the tracker edge … … 431 433 432 434 Float_t L; // track path length 433 Float_t ErrorT; // error on the time measurement 434 435 Float_t ErrorT; // error on the time measurement 436 435 437 Float_t D0; // track signed transverse impact parameter 436 438 Float_t ErrorD0; // signed error on the track signed transverse impact parameter 437 439 438 440 Float_t DZ; // track transverse momentum 439 441 Float_t ErrorDZ; // track transverse momentum error 440 442 441 443 Float_t P; // track transverse momentum 442 444 Float_t ErrorP; // track transverse momentum error 443 445 444 446 Float_t PT; // track transverse momentum 445 447 Float_t ErrorPT; // track transverse momentum error 446 448 447 449 Float_t CtgTheta; // track transverse momentum 448 450 Float_t ErrorCtgTheta; // track transverse momentum error 449 451 450 452 Float_t Phi; // track azimuthal angle 451 453 Float_t ErrorPhi; // track azimuthal angle 452 454 453 455 Float_t Xd; // X coordinate of point of closest approach to vertex 454 456 Float_t Yd; // Y coordinate of point of closest approach to vertex … … 458 460 459 461 Int_t VertexIndex; // reference to vertex 460 462 461 463 static CompBase *fgCompare; //! 462 464 const CompBase *GetCompare() const { return fgCompare; } … … 584 586 585 587 // tracking resolution 586 588 587 589 Float_t TrackResolution; 588 590 … … 612 614 613 615 // vertex variables 614 616 615 617 Int_t ClusterIndex; 616 618 Int_t ClusterNDF; -
modules/PileUpMerger.cc
r0f73bea r5496767 134 134 dz0 = -1.0e6; 135 135 dt0 = -1.0e6; 136 136 137 137 dt *= c_light*1.0E3; // necessary in order to make t in mm/c 138 138 dz *= 1.0E3; // necessary in order to make z in mm 139 139 140 140 //cout<<dz<<","<<dt<<endl; 141 141 142 142 vx = 0.0; 143 143 vy = 0.0; 144 144 145 145 numberOfParticles = fInputArray->GetEntriesFast(); 146 146 nch = 0; 147 sumpt2 = 0.0; 148 147 sumpt2 = 0.0; 148 149 149 while((candidate = static_cast<Candidate*>(fItInputArray->Next()))) 150 150 { … … 154 154 t = candidate->Position.T(); 155 155 pt = candidate->Momentum.Pt(); 156 156 157 157 // take postion and time from first stable particle 158 158 if (dz0 < -999999.0) … … 164 164 candidate->Position.SetZ(z - dz0 + dz); 165 165 candidate->Position.SetT(t - dt0 + dt); 166 166 167 167 fParticleOutputArray->Add(candidate); 168 168 169 169 if(TMath::Abs(candidate->Charge) > 1.0E-9) 170 170 { 171 nch++; 171 nch++; 172 172 sumpt2 += pt*pt; 173 } 173 } 174 174 } 175 175 176 176 if(numberOfParticles > 0) 177 177 { … … 179 179 vy /= sumpt2; 180 180 } 181 181 182 182 nvtx++; 183 183 factory = GetFactory(); … … 188 188 vertex->ClusterNDF = nch; 189 189 vertex->SumPT2 = sumpt2; 190 vertex->GenSumPT2 = sumpt2; 191 190 vertex->GenSumPT2 = sumpt2; 191 192 192 fVertexOutputArray->Add(vertex); 193 193 … … 204 204 case 2: 205 205 numberOfEvents = fMeanPileUp; 206 break; 206 break; 207 207 default: 208 208 numberOfEvents = gRandom->Poisson(fMeanPileUp); … … 233 233 vx = 0.0; 234 234 vy = 0.0; 235 235 236 236 numberOfParticles = 0; 237 237 sumpt2 = 0.0; 238 238 239 vertex = factory->NewCandidate(); 240 239 241 while(fReader->ReadParticle(pid, x, y, z, t, px, py, pz, e)) 240 242 { … … 262 264 vx += candidate->Position.X(); 263 265 vy += candidate->Position.Y(); 264 266 265 267 ++numberOfParticles; 266 268 if(TMath::Abs(candidate->Charge) > 1.0E-9) 267 269 { 268 nch++; 270 nch++; 269 271 sumpt2 += pt*pt; 270 } 271 272 vertex->AddCandidate(candidate); 273 274 } 275 272 276 fParticleOutputArray->Add(candidate); 273 277 } … … 278 282 vy /= numberOfParticles; 279 283 } 280 284 281 285 nvtx++; 282 286 283 vertex = factory->NewCandidate();284 287 vertex->Position.SetXYZT(vx, vy, dz, dt); 285 286 vertex->ClusterIndex = nvtx; 288 289 vertex->ClusterIndex = nvtx; 287 290 vertex->ClusterNDF = nch; 288 291 vertex->SumPT2 = sumpt2; 289 292 vertex->GenSumPT2 = sumpt2; 290 293 291 294 vertex->IsPU = 1; 292 295 293 296 fVertexOutputArray->Add(vertex); 294 297 295 298 } 296 299 } -
modules/TreeWriter.cc
r0f73bea r5496767 240 240 { 241 241 TIter iterator(array); 242 Candidate *candidate = 0 ;242 Candidate *candidate = 0, *constituent = 0; 243 243 Vertex *entry = 0; 244 244 245 245 const Double_t c_light = 2.99792458E8; 246 246 247 247 Double_t x, y, z, t, xError, yError, zError, sigma, sumPT2, btvSumPT2, genDeltaZ, genSumPT2; 248 248 UInt_t index, ndf; 249 249 250 250 array->Sort(); 251 251 252 252 // loop over all vertices 253 253 iterator.Reset(); 254 254 while((candidate = static_cast<Candidate*>(iterator.Next()))) 255 255 { 256 256 257 257 index = candidate->ClusterIndex; 258 258 ndf = candidate->ClusterNDF; … … 262 262 genDeltaZ = candidate->GenDeltaZ; 263 263 genSumPT2 = candidate->GenSumPT2; 264 264 265 265 x = candidate->Position.X(); 266 266 y = candidate->Position.Y(); 267 267 z = candidate->Position.Z(); 268 268 t = candidate->Position.T()*1.0E-3/c_light; 269 269 270 270 xError = candidate->PositionError.X (); 271 271 yError = candidate->PositionError.Y (); … … 281 281 entry->GenDeltaZ = genDeltaZ; 282 282 entry->GenSumPT2 = genSumPT2; 283 283 284 284 entry->X = x; 285 285 entry->Y = y; 286 286 entry->Z = z; 287 287 entry->T = t; 288 288 289 289 entry->ErrorX = xError; 290 290 entry->ErrorY = yError; 291 291 entry->ErrorZ = zError; 292 293 TIter itConstituents(candidate->GetCandidates()); 294 295 while((constituent = static_cast<Candidate*>(itConstituents.Next()))) 296 { 297 entry->Constituents.Add(constituent); 298 } 299 292 300 } 293 301 } … … 332 340 entry->ZOuter = position.Z(); 333 341 entry->TOuter = position.T()*1.0E-3/c_light; 334 342 335 343 entry->L = candidate->L; 336 344 337 345 entry->D0 = candidate->D0; 338 346 entry->ErrorD0 = candidate->ErrorD0; … … 346 354 entry->ErrorCtgTheta = candidate->ErrorCtgTheta; 347 355 entry->Phi = candidate->Phi; 348 entry->ErrorPhi = candidate->ErrorPhi; 349 356 entry->ErrorPhi = candidate->ErrorPhi; 357 350 358 entry->Xd = candidate->Xd; 351 359 entry->Yd = candidate->Yd; … … 361 369 362 370 entry->Eta = eta; 363 371 364 372 particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0)); 365 373 const TLorentzVector &initialPosition = particle->Position;
Note:
See TracChangeset
for help on using the changeset viewer.