source: trunk/source/processes/hadronic/models/util/src/G4SampleResonance.cc@ 1347

Last change on this file since 1347 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 4.8 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28// ------------------------------------------------------------
29// GEANT 4 class header file
30//
31// ---------------- G4SampleResonance ----------------
32// by Henning Weber, March 2001.
33// helper class for sampling resonance masses
34// ------------------------------------------------------------
35
36
37#include "globals.hh"
38#include <iostream>
39#include "G4SampleResonance.hh"
40#include "G4DecayTable.hh"
41#include "Randomize.hh"
42#include "G4HadronicException.hh"
43
44G4SampleResonance::minMassMapType G4SampleResonance::minMassCache;
45
46G4double G4SampleResonance::GetMinimumMass(const G4ParticleDefinition* p) const
47{
48
49 G4double minResonanceMass = DBL_MAX;
50
51 if ( p->IsShortLived() )
52 {
53 minMassMapIterator i = minMassCache.find(p);
54 if ( i!=minMassCache.end() )
55 {
56 minResonanceMass = (*i).second;
57 }
58 else
59 {
60 // G4cout << "--> request for " << p->GetParticleName() << G4endl;
61
62 const G4DecayTable* theDecays = const_cast<G4ParticleDefinition*>(p)->GetDecayTable();
63 const G4int nDecays = theDecays->entries();
64
65 for (G4int i=0; i<nDecays; i++)
66 {
67 const G4VDecayChannel* aDecay = theDecays->GetDecayChannel(i);
68 const G4int nDaughters = aDecay->GetNumberOfDaughters();
69
70 G4double minChannelMass = 0;
71
72 for (G4int j=0; j<nDaughters; j++)
73 {
74 const G4ParticleDefinition* aDaughter = const_cast<G4VDecayChannel*>(aDecay)->GetDaughter(j);
75 G4double minMass = GetMinimumMass(aDaughter);
76 if (!minMass) minMass = DBL_MAX; // exclude gamma channel;
77 minChannelMass+=minMass;
78 }
79 // G4cout << "channel mass for the above is " << minChannelMass/MeV << G4endl;
80 if (minChannelMass < minResonanceMass) minResonanceMass = minChannelMass;
81
82 }
83 // replace this as soon as the compiler supports mutable!!
84 G4SampleResonance* self = const_cast<G4SampleResonance*>(this);
85 (self->minMassCache)[p] = minResonanceMass;
86
87 }
88 }
89 else
90 {
91
92 minResonanceMass = p->GetPDGMass();
93
94 }
95 // G4cout << "minimal mass for " << p->GetParticleName() << " is " << minResonanceMass/MeV << G4endl;
96
97 return minResonanceMass;
98}
99
100
101
102G4double G4SampleResonance::SampleMass(const G4ParticleDefinition* p, const G4double maxMass) const
103{
104 return SampleMass(p->GetPDGMass(), p->GetPDGWidth(), GetMinimumMass(p), maxMass);
105}
106
107
108G4double G4SampleResonance::SampleMass(const G4double poleMass,
109 const G4double gamma,
110 const G4double minMass,
111 const G4double maxMass) const
112{
113 // Chooses a mass randomly between minMass and maxMass
114 // according to a Breit-Wigner function with constant
115 // width gamma and pole poleMass
116
117 if ( minMass > maxMass )
118 {
119 throw G4HadronicException(__FILE__, __LINE__,
120 "SampleResonanceMass: mass range negative (minMass>maxMass)");
121 }
122
123 G4double returnMass;
124
125 if ( gamma < DBL_EPSILON )
126 {
127 returnMass = std::max(minMass, std::min(maxMass, poleMass));
128 }
129 else
130 {
131 double fmin = BrWigInt0(minMass, gamma, poleMass);
132 double fmax = BrWigInt0(maxMass, gamma, poleMass);
133 double f = fmin + (fmax-fmin)*G4UniformRand();
134 returnMass = BrWigInv(f, gamma, poleMass);
135 }
136
137 return returnMass;
138}
139
140
141
142
Note: See TracBrowser for help on using the repository browser.