Changeset acd0621 in git
- Timestamp:
- May 4, 2016, 5:29:47 PM (9 years ago)
- Branches:
- ImprovedOutputFile, Timing, dual_readout, llp, master
- Children:
- 0e2f49b
- Parents:
- 76c2a3b
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
classes/DelphesClasses.h
r76c2a3b racd0621 147 147 Float_t Pz; // particle momentum vector (z component) | hepevt.phep[number][2] 148 148 149 Float_t PT; // particle transverse momentum 149 Float_t D0; 150 Float_t DZ; 151 Float_t P; 152 Float_t PT; 153 Float_t CtgTheta; 154 Float_t Phi; 150 155 Float_t Eta; // particle pseudorapidity 151 Float_t Phi; // particle azimuthal angle152 153 156 Float_t Rapidity; // particle rapidity 154 157 … … 397 400 Int_t Charge; // track charge 398 401 399 Float_t PT; // track transverse momentum400 401 402 Float_t Eta; // track pseudorapidity 402 Float_t Phi; // track azimuthal angle 403 403 404 404 Float_t EtaOuter; // track pseudorapidity at the tracker edge 405 405 Float_t PhiOuter; // track azimuthal angle at the tracker edge … … 415 415 Float_t TOuter; // track position (z component) at the tracker edge 416 416 417 Float_t L; // track path length 418 417 419 Float_t D0; // track signed transverse impact parameter 418 420 Float_t ErrorD0; // signed error on the track signed transverse impact parameter 421 422 Float_t DZ; // track transverse momentum 423 Float_t ErrorDZ; // track transverse momentum error 424 425 Float_t P; // track transverse momentum 426 Float_t ErrorP; // track transverse momentum error 427 428 Float_t PT; // track transverse momentum 429 Float_t ErrorPT; // track transverse momentum error 430 431 Float_t CtgTheta; // track transverse momentum 432 Float_t ErrorCtgTheta; // track transverse momentum error 433 434 Float_t Phi; // track azimuthal angle 435 Float_t ErrorPhi; // track azimuthal angle 436 419 437 Float_t Xd; // X coordinate of point of closest approach to vertex 420 438 Float_t Yd; // Y coordinate of point of closest approach to vertex -
modules/ParticlePropagator.cc
r76c2a3b racd0621 283 283 p = candidateMomentum.P(); 284 284 ctgTheta = 1.0 / TMath::Tan (candidateMomentum.Theta ()); 285 285 286 286 // 3. time evaluation t = TMath::Min(t_r, t_z) 287 287 // t_r : time to exit from the sides … … 337 337 if(r_t > 0.0) 338 338 { 339 340 // store these variables before cloning 341 342 candidate->D0 = d0*1.0E3; 343 candidate->DZ = dz*1.0E3; 344 candidate->P = p; 345 candidate->PT = pt; 346 candidate->CtgTheta = ctgTheta; 347 candidate->Phi = phip; 348 339 349 mother = candidate; 340 350 candidate = static_cast<Candidate*>(candidate->Clone()); … … 344 354 candidate->Momentum = candidateMomentum; 345 355 346 candidate->L = l*1.0E3; 347 candidate->D0 = d0*1.0E3; 348 candidate->DZ = dz*1.0E3; 349 candidate->P = p; 350 candidate->PT = pt; 351 candidate->CtgTheta = ctgTheta; 352 candidate->Phi = phi; 353 354 candidate->Xd = xd*1.0E3; 356 candidate->L = l*1.0E3; 357 358 candidate->Xd = xd*1.0E3; 355 359 candidate->Yd = yd*1.0E3; 356 360 candidate->Zd = zd*1.0E3; -
modules/TreeWriter.cc
r76c2a3b racd0621 215 215 entry->Pz = momentum.Pz(); 216 216 217 entry->D0 = candidate->D0; 218 entry->DZ = candidate->DZ; 219 entry->P = candidate->P; 220 entry->PT = candidate->PT; 221 entry->CtgTheta = candidate->CtgTheta; 222 entry->Phi = candidate->Phi; 223 217 224 entry->Eta = eta; 218 225 entry->Phi = momentum.Phi(); … … 291 298 entry->ZOuter = position.Z(); 292 299 entry->TOuter = position.T()*1.0E-3/c_light; 293 294 entry->D0 = candidate->D0; 295 entry->ErrorD0 = candidate->ErrorD0 ; 300 301 entry->L = candidate->L; 302 303 entry->D0 = candidate->D0; 304 entry->ErrorD0 = candidate->ErrorD0; 305 entry->DZ = candidate->DZ; 306 entry->ErrorDZ = candidate->ErrorDZ; 307 entry->P = candidate->P; 308 entry->ErrorP = candidate->ErrorP; 309 entry->PT = candidate->PT; 310 entry->ErrorPT = candidate->ErrorPT; 311 entry->CtgTheta = candidate->CtgTheta; 312 entry->ErrorCtgTheta = candidate->ErrorCtgTheta; 313 entry->Phi = candidate->Phi; 314 entry->ErrorPhi = candidate->ErrorPhi; 315 296 316 entry->Xd = candidate->Xd; 297 317 entry->Yd = candidate->Yd; … … 307 327 308 328 entry->Eta = eta; 309 entry->Phi = momentum.Phi(); 310 entry->PT = pt; 311 329 312 330 particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0)); 313 331 const TLorentzVector &initialPosition = particle->Position;
Note:
See TracChangeset
for help on using the changeset viewer.