MyEvent.cxx

Go to the documentation of this file.
00001 // Author: Bertrand Bellenot   22/08/02
00002 
00003 /*************************************************************************
00004  * Copyright (C) 1995-2002, Bertrand Bellenot.                           *
00005  * All rights reserved.                                                  *
00006  *                                                                       *
00007  * For the licensing terms see the LICENSE file.                         *
00008  *************************************************************************/
00009 
00010 #include <stdlib.h>
00011 
00012 #include <TROOT.h>
00013 #include <TPolyLine3D.h>
00014 #include <TRandom.h>
00015 #include <TParticle.h>
00016 #include <TDecayChannel.h>
00017 #include <TGenPhaseSpace.h>
00018 #include "MyParticle.h"
00019 #include "TVector3.h"
00020 #include "MyEvent.h"
00021 #include "RootShower.h"
00022 
00023 //______________________________________________________________________________
00024 //
00025 // MyEvent class implementation
00026 //______________________________________________________________________________
00027 
00028 ClassImp(EventHeader)
00029 ClassImp(MyEvent)
00030 
00031 TClonesArray *MyEvent::fgParticles = 0;
00032 
00033 //______________________________________________________________________________
00034 MyEvent::MyEvent()
00035 {
00036    // Create an Event object.
00037    // When the constructor is invoked for the first time, the
00038    // class static variables fgParticles and fgTracks is 0 and
00039    // the TClonesArray fgParticles is created.
00040 
00041    if (!fgParticles) fgParticles = new TClonesArray("MyParticle", 1000);
00042    fgParticles->BypassStreamer(kFALSE);
00043    fParticles = fgParticles;
00044    fNparticles = 0;
00045 }
00046 
00047 //______________________________________________________________________________
00048 MyEvent::~MyEvent()
00049 {
00050    // Destructor
00051 
00052    Clear();
00053 }
00054 
00055 //______________________________________________________________________________
00056 void MyEvent::Init(Int_t id, Int_t first_particle, Double_t E_0, Double_t B_0)
00057 {
00058    // Initialize event ...
00059    // creates detector and set initial values
00060 
00061    Char_t  strtmp[80];
00062    Int_t i;
00063    fId = id;
00064    fB = B_0;
00065 
00066    // generate array of energies threshold used
00067    // to give a track color related to the particle
00068    // energy
00069    for (i=0;i<16;i++)
00070       fEThreshold[i] = EMass + E_0 / (2 << i);
00071 
00072    Clear();
00073    Reset();
00074 
00075    if (!fgParticles) fgParticles = new TClonesArray("MyParticle", 1000);
00076    fParticles = fgParticles;
00077    fNparticles = 0;
00078 
00079    fTotalParticles = 0;
00080    fLast = 0;
00081    fAliveParticles = 1;
00082    fMatter = 0;
00083 
00084    fDetector.ClearELoss();
00085 
00086    TVector3 location(0.0,fDetector.GetMinY(),0.0);
00087    TVector3 momentum(0.0,E_0,0.0);
00088 
00089    AddParticle(0,first_particle, location, momentum);
00090    GetParticle(0)->GenerateTimeOfDecay();
00091 
00092    gTmpLTI = gEventListTree->AddItem(gBaseLTI, GetParticle(0)->GetName());
00093    gTmpLTI->SetUserData(GetParticle(0));
00094    sprintf(strtmp,"%1.2e GeV",GetParticle(0)->Energy());
00095    gEventListTree->SetToolTipItem(gTmpLTI, strtmp);
00096    gLTI[0] = gTmpLTI;
00097 
00098 }
00099 
00100 //______________________________________________________________________________
00101 void MyEvent::Clear(Option_t *option)
00102 {
00103    // Clear tracks and particles arrays
00104 
00105    fgParticles->Delete(option);
00106    fNparticles = 0;
00107    fMatter = 0;
00108 }
00109 
00110 //______________________________________________________________________________
00111 void MyEvent::Reset(Option_t *option)
00112 {
00113    // Static function to reset all static objects for this event
00114 
00115    fgParticles->Delete(option);
00116    fNparticles = 0;
00117    fMatter = 0;
00118 }
00119 
00120 //______________________________________________________________________________
00121 void MyEvent::SetHeader(Int_t i, Int_t run, TDatime date, Int_t primary,
00122                         Double_t energy)
00123 {
00124    // set event header with event identification and startup parameters
00125    
00126    fEvtHdr.Set(i, run, date, primary, energy);
00127 }
00128 
00129 //______________________________________________________________________________
00130 MyParticle *MyEvent::AddParticle(Int_t id, Int_t pdg_code, const TVector3 &pos,
00131                                  const TVector3 &mom)
00132 {
00133    // Add a new particle to the list of particles for this event.
00134    // To avoid calling the very time consuming operator new for each track,
00135    // the standard but not well know C++ operator "new with placement"
00136    // is called. If particle[i] is 0, a new MyParticle object will be created
00137    // otherwise the previous particle[i] will be overwritten.
00138    
00139    TClonesArray &parts = *fParticles;
00140    MyParticle *part = new(parts[fNparticles++])
00141                           MyParticle(id, pdg_code, CREATED, UNDEFINE, pos, mom);
00142    part->AddTrack(pos.x(), pos.y(), pos.z(), ParticleColor(id));
00143    //Save reference to last particle in the collection of particles
00144    fLastParticle = part;
00145    return part;
00146 }
00147 
00148 //______________________________________________________________________________
00149 Int_t MyEvent::Action(Int_t id)
00150 {
00151    // main event's action
00152    
00153    Int_t  nchild;
00154    CheckMatter(id);
00155    if (GetParticle(id)->GetDecayType() == UNDEFINE)
00156       DefineDecay(id);
00157    if (GetParticle(id)->GetPdgCode() == PHOTON) {
00158       // compute the step delta x to be covered by the particle
00159       TVector3 delta_x(GetParticle(id)->GetvMoment() *
00160                (CSpeed * fDetector.GetdT(fMatter) / GetParticle(id)->Energy()));
00161       // check if moved too far (out of detector's bouds)
00162       if (Move(id, delta_x) == DEAD)
00163          // set its status as dead
00164          DeleteParticle(id);
00165       else {
00166          // if distance covered is greater than particle's decay length,
00167          // apply pair production and check if particle is dead. If not,
00168          // increment total alive particles by the two created children,
00169          // then set the particle status as dead
00170          if (GetParticle(id)->GetPassed() >= GetParticle(id)->GetDecayLength()) {
00171             if (PairCreation(id) == DEAD) return(DEAD);
00172             else {
00173                fAliveParticles += 2;
00174                DeleteParticle(id);
00175             }
00176          }
00177       }
00178    }
00179    else if ((GetParticle(id)->GetPdgCode() == NEUTRINO_E) ||
00180            (GetParticle(id)->GetPdgCode() == NEUTRINO_TAU) ||
00181            (GetParticle(id)->GetPdgCode() == NEUTRINO_MUON) ||
00182            (GetParticle(id)->GetPdgCode() == ANTINEUTRINO_E) ||
00183            (GetParticle(id)->GetPdgCode() == ANTINEUTRINO_TAU) ||
00184            (GetParticle(id)->GetPdgCode() == ANTINEUTRINO_MUON) ) {
00185       // if current particle is a neutrino ( or antineutrino )
00186       // set its status as dead ( estimate its probability of
00187       // interaction as null )
00188       DeleteParticle(id);
00189    }
00190    else { // particle is not a photon or neutrino
00191       // if current particle is charged, apply magnetic field influence
00192       if (GetParticle(id)->GetPDG()->Charge() != 0)
00193          ScatterAngle(id);
00194       if ((fB != 0) && (GetParticle(id)->GetPDG()->Charge() != 0))
00195          MagneticField(id);
00196       // compute the step delta x to be covered by the particle
00197       TVector3 delta_x(GetParticle(id)->GetvMoment() * (CSpeed *
00198                        fDetector.GetdT(fMatter) / GetParticle(id)->Energy()));
00199       // check if moved too far (out of detector's bouds)
00200       if (Move(id, delta_x) == DEAD) {
00201          // set its status as dead
00202          DeleteParticle(id);
00203       }
00204       else {
00205          // check energy loss, and if too much energy loss
00206          // ( particle at rest ) set its status as dead
00207          if (DEDX(id) == DEAD) {
00208             DeleteParticle(id);
00209          }
00210          else {
00211             // if at end of particle's life time, decay it
00212             if (CheckDecayTime(id) == 1) {
00213                // if no child found
00214                if ((nchild = Decay(id)) == -1) {
00215                   return(DEAD);
00216                }
00217                else {
00218                   // else increment total alive particles by amount
00219                   // of particle's children
00220                   fAliveParticles += nchild;
00221                   DeleteParticle(id);
00222                }
00223             }
00224             // if not at end of particle's life time, check if distance
00225             // covered is greater than particle's decay length, apply
00226             // defined decay type and check if particle is dead. If not,
00227             // increment total alive particles by the two created children,
00228             // then set the particle status as dead
00229             else if (GetParticle(id)->GetPassed() >=
00230                      GetParticle(id)->GetDecayLength()) {
00231                switch (GetParticle(id)->GetDecayType()) {
00232                   case BREMS:
00233                      if (Bremsstrahlung(id) == DEAD) return(DEAD);
00234                      else {
00235                         fAliveParticles += 2;
00236                         DeleteParticle(id);
00237                      }
00238                      break;
00239                   case CONVERSION:
00240                      if (PairCreation(id) == DEAD) return(DEAD);
00241                      else {
00242                         fAliveParticles += 2;
00243                         DeleteParticle(id);
00244                      }
00245                      break;
00246                }
00247             }
00248          }
00249       }
00250    }
00251    return(ALIVE);
00252 }
00253 
00254 //______________________________________________________________________________
00255 Double_t MyEvent::BremsProb(Int_t id)
00256 {
00257    // Check if bremsstrahlung is allowed and generate
00258    // a random decay length related to detector's material
00259    // radiation length (X0)
00260    
00261    Double_t p, retval;
00262 
00263    if (GetParticle(id)->Energy() > GetParticle(id)->GetMass()) {
00264       p = gRandom->Uniform(0.0, 1.0);
00265       retval = (-fDetector.GetX0(fMatter))*TMath::Log(p);
00266       return (retval);
00267    }
00268    else return (-1.);
00269 }
00270 
00271 //______________________________________________________________________________
00272 Int_t MyEvent::Bremsstrahlung(Int_t id)
00273 {
00274     // compute bremsstrahlung for particle "id"
00275     
00276    Double_t  ratio;
00277    Int_t     d_num1,d_num2;
00278    Char_t    strtmp[80];
00279    MyParticle *part;
00280 
00281    // find two ids for children particles
00282    if ((FindFreeId(&d_num1) != DEAD) && (FindFreeId(&d_num2) != DEAD)) {
00283       // compute the particle's energy ratio...
00284       ratio = (GetParticle(id)->Energy() - GetParticle(id)->GetMass()) /
00285               (2 * GetParticle(id)->P());
00286       // create first child if fact, electron continues with less energy
00287       // and in a different direction. To that end the electron is added
00288       // to its own list of children, because otherwise it would vanish.
00289       part = AddParticle(d_num1, GetParticle(id)->GetPdgCode(),
00290                          GetParticle(id)->GetvLocation(),
00291                          GetParticle(id)->GetvMoment() * ratio);
00292       part->SetFirstMother(id);
00293       // as its first child is in fact the same particle,
00294       // keep the same decay time
00295       part->SetTimeOfDecay(GetParticle(id)->GetTimeOfDecay());
00296       GetParticle(id)->SetChild(0, d_num1);
00297       // add a track related to this particle
00298 
00299       // add a particle related list tree item to the event list tree
00300       gTmpLTI = gEventListTree->AddItem(gLTI[id], part->GetName());
00301       gTmpLTI->SetUserData(part);
00302       sprintf(strtmp,"%1.2e GeV",part->Energy());
00303       gEventListTree->SetToolTipItem(gTmpLTI, strtmp);
00304       gLTI[d_num1] = gTmpLTI;
00305 
00306       // create second child
00307       part = AddParticle(d_num2,PHOTON, GetParticle(id)->GetvLocation(),
00308                          GetParticle(id)->GetvMoment() * ratio);
00309       part->SetFirstMother(id);
00310       // generate time of decay (not used in this case, as it is a photon,
00311       // but to keep the same philosophy in every case...
00312       part->GenerateTimeOfDecay();
00313       GetParticle(id)->SetChild(1, d_num2);
00314 
00315       // add a track related to this particle
00316 
00317       // add a particle related list tree item to the event list tree
00318       gTmpLTI = gEventListTree->AddItem(gLTI[id], part->GetName());
00319       gTmpLTI->SetUserData(part);
00320       sprintf(strtmp,"%1.2e GeV",part->Energy());
00321       gEventListTree->SetToolTipItem(gTmpLTI, strtmp);
00322       gLTI[d_num2] = gTmpLTI;
00323 
00324       // increment number of children by the two created particles
00325       GetParticle(id)->SetNChildren(2);
00326 
00327       return(ALIVE);
00328    }
00329    else return(DEAD);
00330 }
00331 
00332 //______________________________________________________________________________
00333 Int_t MyEvent::CheckDecayTime(Int_t id)
00334 {
00335    // Check decay time for particle "id".
00336 
00337    if ( (GetParticle(id)->GetPdgCode() == PHOTON) ||
00338         (GetParticle(id)->GetPdgCode() == ELECTRON) ||
00339         (GetParticle(id)->GetPdgCode() == POSITRON ))
00340       return 0;
00341    Double_t timeofdecay = GetParticle(id)->GetTimeOfDecay();
00342    if (timeofdecay == 0.0)  return 0;
00343    Double_t distToDecay = timeofdecay * 0.996 * CSpeed;
00344    // check if actual particle life is greater than particle life time
00345    if (GetParticle(id)->GetPassed() >= distToDecay) {
00346       return 1;
00347    }
00348    return 0;
00349 }
00350 
00351 //______________________________________________________________________________
00352 void MyEvent::CheckMatter(Int_t id)
00353 {
00354    // Check material into which the particle "id" is. 
00355 
00356    TGeoNode *Node = gGeoManager->FindNode(
00357              GetParticle(id)->GetvLocation().x(),
00358              GetParticle(id)->GetvLocation().y(),
00359              GetParticle(id)->GetvLocation().z());
00360    if (Node) fMatter = Node->GetNumber();
00361 }
00362 
00363 //______________________________________________________________________________
00364 Int_t MyEvent::Decay(Int_t id)
00365 {
00366    // Decay the particle "id".
00367 
00368    Char_t   strtmp[80];
00369    Int_t    d_num[5];
00370    Int_t    n_daughters;
00371    Int_t    ptype[5];
00372    Double_t mass[5];
00373    Int_t    i, index;
00374    Double_t sumBR = 0.0;
00375    MyParticle *Particle[5];
00376    MyParticle *part;
00377 
00378    // compute total branching ratio
00379    for (i=0;i<GetParticle(id)->GetPDG()->NDecayChannels();i++) {
00380       sumBR += GetParticle(id)->GetPDG()->DecayChannel(i)->BranchingRatio();
00381    }
00382    // choose random decay in respect to the branching ratio
00383 again:
00384    float r = gRandom->Uniform(sumBR);
00385    index = 0;
00386    while
00387       ((r -= GetParticle(id)->GetPDG()->DecayChannel(index)->BranchingRatio())
00388        > 0 && index < GetParticle(id)->GetPDG()->NDecayChannels()) index++;
00389 
00390    // set number of daughters
00391    n_daughters = GetParticle(id)->GetPDG()->DecayChannel(index)->NDaughters();
00392    for (i=0;i<n_daughters;i++) {
00393       // create temporary child particle to obtain its mass
00394       ptype[i] =
00395          GetParticle(id)->GetPDG()->DecayChannel(index)->DaughterPdgCode(i);
00396       if (TMath::Abs(ptype[i]) < 6) // it is a quark...do it again
00397          goto again;
00398       Particle[i] = new MyParticle(0,ptype[i], CREATED, UNDEFINE,
00399                                    GetParticle(id)->GetvLocation(),
00400                                    GetParticle(id)->GetvMoment());
00401       mass[i] = Particle[i]->GetMass();
00402       delete Particle[i];
00403    }
00404 
00405    // setup the decay
00406    TLorentzVector W(GetParticle(id)->GetvMoment(), GetParticle(id)->Energy());
00407    TGenPhaseSpace genPhaseSpace;
00408    if (!genPhaseSpace.SetDecay( W, n_daughters, mass ))
00409       return (-1);
00410    genPhaseSpace.Generate();
00411 
00412    // find ids for children
00413    for (i=0;i<n_daughters;i++) {
00414       if (FindFreeId(&d_num[i]) == DEAD) return -1;
00415    }
00416 
00417    TLorentzVector *p;
00418    for (i=0;i<n_daughters;i++) {
00419       p = genPhaseSpace.GetDecay(i);
00420       // create child
00421       part = AddParticle(d_num[i], ptype[i], GetParticle(id)->GetvLocation(),
00422                          p->Vect());
00423       part->SetFirstMother(id);
00424       // generate time of decay (may be useful in this case)
00425       part->GenerateTimeOfDecay();
00426       GetParticle(id)->SetChild(i, d_num[i]);
00427       // add a track related to this child
00428 
00429       // add a child related list tree item to the event list tree
00430       gTmpLTI = gEventListTree->AddItem(gLTI[id], part->GetName());
00431       gTmpLTI->SetUserData(part);
00432       sprintf(strtmp,"%1.2e GeV",part->Energy());
00433       gEventListTree->SetToolTipItem(gTmpLTI, strtmp);
00434       gLTI[d_num[i]] = gTmpLTI;
00435    }
00436    // increment number of children by the number of created particles
00437    GetParticle(id)->SetNChildren(n_daughters);
00438    return(n_daughters);
00439 }
00440 
00441 //______________________________________________________________________________
00442 void MyEvent::DefineDecay(Int_t id)
00443 {
00444    // Define decay type for particle "id", then check decay length for it
00445 
00446    Double_t idecay_length = -1.;
00447    Double_t iactual_length;
00448    Int_t    idecay_type = CONVERSION;
00449 
00450    if ( (GetParticle(id)->GetPdgCode() == ELECTRON) ||
00451         (GetParticle(id)->GetPdgCode() == POSITRON)) {
00452       // check if bremsstrahlung is allowed
00453       if ( (iactual_length = BremsProb(id)) > 0.) {
00454          if ( (idecay_length == -1) || (iactual_length < idecay_length) ) {
00455             idecay_length = iactual_length;
00456             idecay_type = BREMS;
00457          }
00458       }
00459    }
00460    else if (GetParticle(id)->GetPdgCode() == PHOTON) {
00461       // check if pair production is allowed
00462       if ( (iactual_length = PairProb(id)) > 0. ) {
00463          if ( (idecay_length == -1) ||
00464              (iactual_length < idecay_length) ) {
00465             idecay_length = iactual_length;
00466             idecay_type = CONVERSION;
00467          }
00468       }
00469    }
00470    if ( idecay_length > 0) {
00471       GetParticle(id)->SetDecayType(idecay_type);
00472       GetParticle(id)->SetDecayLength(idecay_length);
00473    }
00474    else {
00475       GetParticle(id)->SetDecayType(STABLE);
00476       GetParticle(id)->SetDecayLength(0.0);
00477    }
00478 }
00479 
00480 //______________________________________________________________________________
00481 void MyEvent::DeleteParticle(Int_t id)
00482 {
00483    // Delete the particle id
00484 
00485    // To be sure that the last track has at least two points
00486    if (GetParticle(id)->GetTrack(GetParticle(id)->GetNTracks()-1)->GetN() < 2) {
00487       GetParticle(id)->SetNextPoint(GetParticle(id)->GetTrack(GetParticle(id)->GetNTracks()-1)->GetLineColor());
00488    }
00489    // Add this particle's energy loss at the total
00490    // energy loss into the detector
00491    fDetector.AddELoss(GetParticle(id)->GetELoss());
00492    // Mark the particle's status as dead and decrement
00493    // the total alive particles
00494    GetParticle(id)->SetStatus(DEAD);
00495    fAliveParticles --;
00496 }
00497 
00498 //______________________________________________________________________________
00499 Int_t MyEvent::DEDX(Int_t id)
00500 {
00501    // Compute de/dx for particle "id" into detector material
00502    // for more infos, please refer to the particle data booklet
00503    // from which the formulas has been extracted
00504 
00505    Double_t gamma,abs_beta,abs_p,abs_loss,dX;
00506 
00507    // if particle's energy is equal to its mass, it is at rest,
00508    // so set its status as dead
00509    if (GetParticle(id)->Energy() <= GetParticle(id)->GetMass()) {
00510       return(DEAD);
00511    }
00512    else {
00513       // absolute value of momentum
00514       abs_p = GetParticle(id)->P();
00515       if (abs_p <= 0) {
00516          // if absolute value of momentum is less or equal to zero,
00517          // set it to the particle's mass (minimum allowed value for momentum)
00518          GetParticle(id)->SetMomentum(0.0, 0.0, 0.0,GetParticle(id)->GetMass());
00519       }
00520       else {
00521          // Compute energy loss in detector's material
00522          // cf Bethe Bloch formula
00523          TVector3 p_0(GetParticle(id)->GetvMoment() * (1 / abs_p));
00524          abs_beta = abs_p / GetParticle(id)->Energy();
00525          dX = fDetector.GetdT(fMatter) * CSpeed * abs_beta;
00526          abs_beta *= abs_beta;
00527          if (abs_beta < .9999999999) gamma = 1/TMath::Sqrt(1.0-abs_beta);
00528          else gamma = MAX_GAMMA;
00529          abs_loss = (fDetector.GetPreconst(fMatter) * dX / abs_beta) *
00530                     (TMath::Log(2.0 * GetParticle(id)->GetMass() *
00531                      gamma * gamma * abs_beta /
00532                      fDetector.GetI(fMatter)) - abs_beta);
00533          if (abs_loss < 0) abs_loss = -abs_loss;
00534          if (abs_loss >= (GetParticle(id)->Energy() -
00535                          GetParticle(id)->GetMass())) {
00536             // if energy loss leave less energy to the particle than
00537             // its mass, set its momentum equal to its mass
00538             // (minimum allowed value for momentum)
00539             GetParticle(id)->SetMomentum(0.0, 0.0, 0.0,
00540                              GetParticle(id)->GetMass());
00541          }
00542          else {
00543             // else decrease its energy by calculated energy loss
00544             GetParticle(id)->SetMoment(GetParticle(id)->GetvMoment(),
00545                 GetParticle(id)->Energy() - abs_loss);
00546             abs_p = TMath::Sqrt((GetParticle(id)->Energy() *
00547                           GetParticle(id)->Energy()) -
00548                          (GetParticle(id)->GetMass() *
00549                           GetParticle(id)->GetMass()));
00550             GetParticle(id)->SetMoment(p_0 * abs_p);
00551             // Add calculated energy loss at total particle's energy loss
00552             GetParticle(id)->AddELoss(abs_loss);
00553          }
00554       }
00555       if (GetParticle(id)->Energy() > GetParticle(id)->GetMass()) {
00556          return(ALIVE);
00557       }
00558       else {
00559          return(DEAD);
00560       }
00561    }
00562 }
00563 
00564 //______________________________________________________________________________
00565 Int_t MyEvent::FindFreeId(Int_t *FreeId)
00566 {
00567    // Give next available particle's id.
00568    
00569    fTotalParticles++;
00570    *FreeId = fTotalParticles;
00571    if (fTotalParticles > fLast) fLast = fTotalParticles;
00572    return(ALIVE);
00573 }
00574 
00575 //______________________________________________________________________________
00576 void MyEvent::MagneticField(Int_t id)
00577 {
00578    // Extrapolate track in a constant field oriented along X axis
00579    // translated to C++ from GEANT3 routine GHELX3.
00580    
00581    Double_t sint, sintt, tsint, cos1t, sin2;
00582    Double_t f1, f2, f3, v1, v2, v3;
00583    Double_t pol = GetParticle(id)->GetPDG()->Charge() / 3.0;
00584    Double_t h4  = pol * 2.9979251e-04 * fB;
00585    Double_t tet = -h4 * CSpeed * fDetector.GetdT(fMatter) /
00586                   GetParticle(id)->P();
00587    if (TMath::Abs(tet) > 0.15) {
00588       sint  = TMath::Sin(tet);
00589       sintt = sint / tet;
00590       tsint = (tet - sint) / tet;
00591       sin2  = TMath::Sin(0.5 * tet);
00592       cos1t = 2.0 * sin2 * sin2 / tet;
00593    } else {
00594       tsint = tet * tet / 6.0;
00595       sintt = 1.0 - tsint;
00596       sint  = tet * sintt;
00597       cos1t = 0.5 * tet;
00598    }
00599    f1 = -tet * cos1t;
00600    f2 = sint;
00601    f3 = tet * cos1t * GetParticle(id)->Px();
00602    v1 = GetParticle(id)->Px() + (f1 * GetParticle(id)->Px() + f3);
00603    v2 = GetParticle(id)->Py() + (f1 * GetParticle(id)->Py() + f2 *
00604         GetParticle(id)->Pz());
00605    v3 = GetParticle(id)->Pz() + (f1 * GetParticle(id)->Pz() - f2 *
00606         GetParticle(id)->Py());
00607    TVector3 new_mom(v1, v2, v3);
00608    GetParticle(id)->SetMoment(new_mom);
00609 }
00610 
00611 //______________________________________________________________________________
00612 Int_t MyEvent::Move(Int_t id, TVector3 &dist)
00613 {
00614    // Move particle "id" by step dist, update the distance covered
00615    // then check if out of detector's bounds.
00616 
00617    GetParticle(id)->SetLocation(GetParticle(id)->GetvLocation() + dist);
00618    GetParticle(id)->SetPassed(GetParticle(id)->GetPassed() + dist.Mag());
00619 
00620    if ((GetParticle(id)->GetvLocation().x() > fDetector.GetMaxX()) ||
00621       (GetParticle(id)->GetvLocation().x() < fDetector.GetMinX()) ||
00622       (GetParticle(id)->GetvLocation().y() > fDetector.GetMaxY()) ||
00623       (GetParticle(id)->GetvLocation().y() < fDetector.GetMinY()) ||
00624       (GetParticle(id)->GetvLocation().z() > fDetector.GetMaxZ()) ||
00625       (GetParticle(id)->GetvLocation().z() < fDetector.GetMinZ())) {
00626       return(DEAD);
00627    }
00628    // If not out of bounds, set related Track's next point
00629    else {
00630       if ((GetParticle(id)->GetPdgCode() != PHOTON) &&
00631          (GetParticle(id)->GetPdgCode() != NEUTRINO_E) &&
00632          (GetParticle(id)->GetPdgCode() != NEUTRINO_MUON) &&
00633          (GetParticle(id)->GetPdgCode() != NEUTRINO_TAU) &&
00634          (GetParticle(id)->GetPdgCode() != ANTINEUTRINO_E) &&
00635          (GetParticle(id)->GetPdgCode() != ANTINEUTRINO_MUON) &&
00636          (GetParticle(id)->GetPdgCode() != ANTINEUTRINO_TAU) ) {
00637          GetParticle(id)->SetNextPoint(ParticleColor(id));
00638       }
00639       return(ALIVE);
00640    }
00641 }
00642 
00643 //______________________________________________________________________________
00644 Double_t MyEvent::PairProb(Int_t id)
00645 {
00646    // Check if pair production is allowed and generate
00647    // a random decay length related to detector's material
00648    // radiation length (X0).
00649 
00650    Double_t p;
00651 
00652    if (GetParticle(id)->Energy() > 2.0 * EMass) {
00653       p = gRandom->Uniform(0.0, 1.0);
00654       return ((-9.)*fDetector.GetX0(fMatter)*TMath::Log(p)/7.);
00655    }
00656    return (-1.);
00657 }
00658 
00659 //______________________________________________________________________________
00660 Int_t MyEvent::PairCreation(Int_t id)
00661 {
00662    // Compute the pair production for particle "id"
00663 
00664    Int_t    d_num1, d_num2;
00665    Char_t   strtmp[80];
00666    MyParticle *part;
00667    TGenPhaseSpace genPhaseSpace;
00668 
00669    // setup the decay
00670    TLorentzVector target(0.0, 0.0, 0.0, 0.00511);
00671    TLorentzVector beam(GetParticle(id)->GetvMoment(),
00672                        GetParticle(id)->Energy());
00673    TLorentzVector W = beam + target;
00674    Double_t masses[2] = { EMass, EMass } ;
00675 
00676    if (!genPhaseSpace.SetDecay( W, 2, masses ))
00677       return (DEAD);
00678    genPhaseSpace.Generate();
00679 
00680    TLorentzVector *p;
00681    // calculate the location of the decay vertex
00682 
00683    // find two ids for children particles
00684    if ((FindFreeId(&d_num1) != DEAD) && (FindFreeId(&d_num2) != DEAD)) {
00685       p = genPhaseSpace.GetDecay(0);
00686       // create child
00687       part = AddParticle(d_num1, POSITRON, GetParticle(id)->GetvLocation(),
00688                          p->Vect());
00689       part->SetFirstMother(id);
00690 
00691       // generate time of decay (not used in this case, as it is an electron
00692       // or a positron, but to keep the same philosophy in every case...
00693       part->GenerateTimeOfDecay();
00694       GetParticle(id)->SetChild(0, d_num1);
00695       // add a track related to this particle
00696 
00697       // add a particle related list tree item to the event list tree
00698       gTmpLTI = gEventListTree->AddItem(gLTI[id], part->GetName());
00699       gTmpLTI->SetUserData(part);
00700       sprintf(strtmp,"%1.2e GeV",part->Energy());
00701       gEventListTree->SetToolTipItem(gTmpLTI, strtmp);
00702       gLTI[d_num1] = gTmpLTI;
00703 
00704       // create second child
00705       p = genPhaseSpace.GetDecay(1);
00706 
00707       // create child
00708       part = AddParticle(d_num2, ELECTRON, GetParticle(id)->GetvLocation(),
00709                          p->Vect());
00710       part->SetFirstMother(id);
00711       // generate time of decay (not used in this case, as it is an electron
00712       // or a positron, but to keep the same philosophy in every case...
00713       part->GenerateTimeOfDecay();
00714       GetParticle(id)->SetChild(1, d_num2);
00715       // add a track related to this particle
00716 
00717       // add a particle related list tree item to the event list tree
00718       gTmpLTI = gEventListTree->AddItem(gLTI[id], part->GetName());
00719       gTmpLTI->SetUserData(part);
00720       sprintf(strtmp,"%1.2e GeV",part->Energy());
00721       gEventListTree->SetToolTipItem(gTmpLTI, strtmp);
00722       gLTI[d_num2] = gTmpLTI;
00723 
00724       // increment number of children by the two created particles
00725       GetParticle(id)->SetNChildren(2);
00726 
00727       return(ALIVE);
00728    }
00729    else return(DEAD);
00730 }
00731 
00732 //______________________________________________________________________________
00733 Int_t MyEvent::ParticleColor(Int_t id)
00734 {
00735    // Return color index related to particle's energy.
00736    //Int_t ctable[11] = {2,50,46,45,44,43,42,41,21,19,5};
00737 
00738    Int_t i;
00739    for (i=0;i<16;i++)
00740       if (GetParticle(id)->Energy() > fEThreshold[i]) break;
00741    if (i > 16) i = 16;
00742    return(gColIndex + i);
00743 }
00744 
00745 //______________________________________________________________________________
00746 void MyEvent::ScatterAngle(Int_t id)
00747 {
00748    // Compute scatter angle into the detector's material
00749    // for the current particle
00750    // for more infos, please refer to the particle data booklet
00751    // from which the formulas has been extracted :
00752    // Multiple scattering through small angles.
00753    
00754    Double_t alpha,beta;
00755    Double_t abs_p,p1,p2,r_2;
00756    Double_t fact1,fact2;
00757 
00758    do {
00759       p1 = gRandom->Uniform(-1.0, 1.0);
00760       p2 = gRandom->Uniform(-1.0, 1.0);
00761       r_2 = (p1 * p1) + (p2 * p2);
00762    } while (r_2 > 1.0);
00763    abs_p = GetParticle(id)->P();
00764    alpha = TMath::Sqrt(-2.0 * TMath::Log(r_2) / r_2) *
00765            fDetector.GetTheta0(fMatter) / abs_p;
00766    beta  = gRandom->Uniform(0.0, 2.0 * TMath::Pi());
00767    alpha *= p1;
00768    TVector3 x_0(GetParticle(id)->GetvMoment().Orthogonal());
00769    TVector3 p_0(GetParticle(id)->GetvMoment() * (1.0 / abs_p));
00770    TVector3 y_0(x_0.Cross(p_0));
00771    fact1 = TMath::Sin(alpha);
00772    fact2 = fact1 * TMath::Cos(beta);
00773    fact1 *= TMath::Sin(beta);
00774    TVector3 vtmp1(x_0 * fact1);
00775    TVector3 vtmp2(y_0 * fact2);
00776    TVector3 vtmp3(vtmp2 + p_0);
00777 
00778    GetParticle(id)->SetMoment(vtmp1 + vtmp3);
00779    GetParticle(id)->SetMoment(GetParticle(id)->GetvMoment() *
00780                              (abs_p/ GetParticle(id)->P()));
00781 }
00782 

Generated on Tue Jul 5 15:15:11 2011 for ROOT_528-00b_version by  doxygen 1.5.1