G4SampleResonance.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 //
00028 // ------------------------------------------------------------
00029 //      GEANT 4 class header file
00030 //
00031 //      ---------------- G4SampleResonance ----------------
00032 //             by Henning Weber, March 2001.
00033 //      helper class for sampling resonance masses
00034 // ------------------------------------------------------------
00035 
00036 
00037 #include "globals.hh"
00038 #include <iostream>
00039 #include "G4SampleResonance.hh"
00040 #include "G4DecayTable.hh"
00041 #include "Randomize.hh"
00042 #include "G4HadronicException.hh"
00043 
00044 G4SampleResonance::minMassMapType G4SampleResonance::minMassCache;
00045 
00046 G4double G4SampleResonance::GetMinimumMass(const G4ParticleDefinition* p) const
00047 {
00048 
00049         G4double minResonanceMass = DBL_MAX;
00050 
00051         if ( p->IsShortLived() )
00052         {
00053                 minMassMapIterator iter = minMassCache.find(p);
00054                 if ( iter!=minMassCache.end() )
00055                 {
00056                         minResonanceMass = (*iter).second;
00057                 }
00058                 else
00059                 {
00060                         // G4cout << "--> request for " << p->GetParticleName() << G4endl;
00061 
00062                         const G4DecayTable* theDecays = const_cast<G4ParticleDefinition*>(p)->GetDecayTable();
00063                         const G4int nDecays = theDecays->entries();
00064 
00065                         for (G4int i=0; i<nDecays; i++)
00066                         {
00067                                 const G4VDecayChannel* aDecay = theDecays->GetDecayChannel(i);
00068                                 const G4int nDaughters = aDecay->GetNumberOfDaughters();
00069 
00070                                 G4double minChannelMass = 0;
00071 
00072                                 for (G4int j=0; j<nDaughters; j++)
00073                                 {
00074                                         const G4ParticleDefinition* aDaughter = const_cast<G4VDecayChannel*>(aDecay)->GetDaughter(j);
00075                                         G4double minMass = GetMinimumMass(aDaughter);
00076                                         if (!minMass) minMass = DBL_MAX; // exclude gamma channel;
00077                                         minChannelMass+=minMass;
00078                                 }
00079                                 // G4cout << "channel mass for the above is " << minChannelMass/MeV << G4endl;
00080                                 if (minChannelMass < minResonanceMass) minResonanceMass = minChannelMass;
00081 
00082                         }
00083                         // replace this as soon as the compiler supports mutable!!
00084                         G4SampleResonance* self = const_cast<G4SampleResonance*>(this);
00085                         (self->minMassCache)[p] = minResonanceMass;
00086 
00087                 }
00088         } 
00089         else
00090         {
00091 
00092                 minResonanceMass = p->GetPDGMass();
00093 
00094         }
00095         // G4cout << "minimal mass for " << p->GetParticleName() << " is " << minResonanceMass/MeV << G4endl;
00096 
00097         return minResonanceMass;
00098 } 
00099 
00100 
00101 
00102 G4double G4SampleResonance::SampleMass(const G4ParticleDefinition* p, const G4double maxMass) const
00103 {
00104         return SampleMass(p->GetPDGMass(), p->GetPDGWidth(), GetMinimumMass(p), maxMass);
00105 }
00106 
00107 
00108 G4double G4SampleResonance::SampleMass(const G4double poleMass, 
00109                 const G4double gamma,
00110                 const G4double minMass,
00111                 const G4double maxMass) const
00112 {
00113         // Chooses a mass randomly between minMass and maxMass
00114         //     according to a Breit-Wigner function with constant
00115         //     width gamma and pole poleMass
00116 
00117         if ( minMass > maxMass )
00118         {
00119                 throw G4HadronicException(__FILE__, __LINE__,
00120                                 "SampleResonanceMass: mass range negative (minMass>maxMass)");
00121         }
00122 
00123         G4double returnMass;
00124 
00125         if ( gamma < DBL_EPSILON )
00126         {
00127                 returnMass = std::max(minMass, std::min(maxMass, poleMass));
00128         }
00129         else
00130         {
00131                 double fmin = BrWigInt0(minMass, gamma, poleMass);
00132                 double fmax = BrWigInt0(maxMass, gamma, poleMass);
00133                 double f = fmin + (fmax-fmin)*G4UniformRand();
00134                 returnMass = BrWigInv(f, gamma, poleMass);
00135         }
00136 
00137         return returnMass;
00138 }
00139 
00140 
00141 
00142 

Generated on Mon May 27 17:49:47 2013 for Geant4 by  doxygen 1.4.7