00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #ifndef G4RadioactiveDecay_h
00027 #define G4RadioactiveDecay_h 1
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00056
00057 #include <vector>
00058 #include <map>
00059 #include <CLHEP/Units/SystemOfUnits.h>
00060
00061 #include "G4ios.hh"
00062 #include "globals.hh"
00063 #include "G4VRestDiscreteProcess.hh"
00064 #include "G4ParticleChangeForRadDecay.hh"
00065
00066
00067 #include "G4NucleusLimits.hh"
00068 #include "G4RadioactiveDecayRate.hh"
00069 #include "G4RadioactiveDecayRateVector.hh"
00070 #include "G4RIsotopeTable.hh"
00071 #include "G4RadioactivityTable.hh"
00072 #include "G4ThreeVector.hh"
00073
00074 class G4RadioactiveDecaymessenger;
00075
00076 typedef std::vector<G4RadioactiveDecayRateVector> G4RadioactiveDecayRateTable;
00077 typedef std::vector<G4RadioactiveDecayRate> G4RadioactiveDecayRates;
00078
00079
00080 class G4RadioactiveDecay : public G4VRestDiscreteProcess
00081 {
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092 public:
00093
00094 G4RadioactiveDecay(const G4String& processName="RadioactiveDecay");
00095 ~G4RadioactiveDecay();
00096
00097 G4bool IsApplicable(const G4ParticleDefinition&);
00098
00099
00100
00101
00102 G4bool IsLoaded(const G4ParticleDefinition &);
00103
00104
00105 void SelectAVolume(const G4String aVolume);
00106
00107
00108 void DeselectAVolume(const G4String aVolume);
00109
00110
00111 void SelectAllVolumes();
00112
00113
00114 void DeselectAllVolumes();
00115
00116
00117 void SetDecayBias (G4String filename);
00118
00119
00120 void SetHLThreshold (G4double hl) {halflifethreshold = hl;}
00121
00122
00123 void SetICM (G4bool icm) {applyICM = icm;}
00124
00125
00126 void SetARM (G4bool arm) {applyARM = arm;}
00127
00128
00129 void SetSourceTimeProfile (G4String filename) ;
00130
00131
00132 G4bool IsRateTableReady(const G4ParticleDefinition &);
00133
00134
00135
00136
00137 void AddDecayRateTable(const G4ParticleDefinition&);
00138
00139
00140
00141
00142
00143 void GetDecayRateTable(const G4ParticleDefinition&);
00144
00145
00146
00147
00148
00149 void SetDecayRate(G4int,G4int,G4double, G4int, std::vector<G4double>,
00150 std::vector<G4double>);
00151
00152
00153
00154 std::vector<G4RadioactivityTable*> GetTheRadioactivityTables()
00155 {return theRadioactivityTables;}
00156
00157
00158 G4DecayTable *LoadDecayTable (G4ParticleDefinition & theParentNucleus);
00159
00160
00161 void AddUserDecayDataFile(G4int Z, G4int A,G4String filename);
00162
00163
00164
00165
00166
00167
00168 inline void SetVerboseLevel(G4int value) {verboseLevel = value;}
00169
00170
00171 inline G4int GetVerboseLevel() const {return verboseLevel;}
00172
00173
00174 inline void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
00175 {theNucleusLimits = theNucleusLimits1 ;}
00176
00177
00178
00179 inline G4NucleusLimits GetNucleusLimits() const
00180 {return theNucleusLimits;}
00181
00182
00183
00184 inline void SetAnalogueMonteCarlo (G4bool r ) {
00185 AnalogueMC = r;
00186 if (!AnalogueMC) halflifethreshold = 1e-6*CLHEP::s;
00187 }
00188
00189
00190
00191 inline void SetFBeta (G4bool r ) { FBeta = r; }
00192
00193
00194 inline G4bool IsAnalogueMonteCarlo () {return AnalogueMC;}
00195
00196
00197
00198 inline void SetBRBias (G4bool r) {
00199 BRBias = r;
00200 SetAnalogueMonteCarlo(0);
00201 }
00202
00203
00204 inline void SetSplitNuclei (G4int r) {
00205 NSplit = r;
00206 SetAnalogueMonteCarlo(0);
00207 }
00208
00209
00210 inline G4int GetSplitNuclei () {return NSplit;}
00211
00212
00213 inline void SetDecayDirection(const G4ThreeVector& theDir) {
00214 forceDecayDirection = theDir.unit();
00215 }
00216
00217 inline const G4ThreeVector& GetDecayDirection() const {
00218 return forceDecayDirection;
00219 }
00220
00221 inline void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg) {
00222 forceDecayHalfAngle = std::min(std::max(0.*CLHEP::deg,halfAngle),180.*CLHEP::deg);
00223 }
00224
00225 inline G4double GetDecayHalfAngle() const {return forceDecayHalfAngle;}
00226
00227 inline void SetDecayCollimation(const G4ThreeVector& theDir,
00228 G4double halfAngle = 0.*CLHEP::deg) {
00229 SetDecayDirection(theDir);
00230 SetDecayHalfAngle(halfAngle);
00231 }
00232
00233
00234
00235
00236 void BuildPhysicsTable(const G4ParticleDefinition &);
00237
00238 protected:
00239
00240 G4VParticleChange* DecayIt(const G4Track& theTrack,
00241 const G4Step& theStep);
00242
00243 G4DecayProducts* DoDecay(G4ParticleDefinition& theParticleDef);
00244
00245
00246 void CollimateDecay(G4DecayProducts* products);
00247 void CollimateDecayProduct(G4DynamicParticle* product);
00248 G4ThreeVector ChooseCollimationDirection() const;
00249
00250 G4double GetMeanFreePath(const G4Track& theTrack, G4double previousStepSize,
00251 G4ForceCondition* condition);
00252
00253 G4double GetMeanLifeTime(const G4Track& theTrack,
00254 G4ForceCondition* condition);
00255
00256 G4double GetTaoTime(const G4double,const G4double);
00257
00258 G4double GetDecayTime();
00259
00260 G4int GetDecayTimeBin(const G4double aDecayTime);
00261
00262 private:
00263
00264 G4RadioactiveDecay(const G4RadioactiveDecay &right);
00265 G4RadioactiveDecay & operator=(const G4RadioactiveDecay &right);
00266
00267 G4RadioactiveDecaymessenger* theRadioactiveDecaymessenger;
00268 G4RIsotopeTable* theIsotopeTable;
00269
00270 G4NucleusLimits theNucleusLimits;
00271
00272 G4double HighestValue;
00273
00274 G4bool isInitialised;
00275 G4bool AnalogueMC;
00276 G4bool BRBias;
00277 G4bool FBeta;
00278 G4int NSplit;
00279
00280 G4double halflifethreshold;
00281 G4bool applyICM;
00282 G4bool applyARM;
00283
00284
00285 G4ThreeVector forceDecayDirection;
00286 G4double forceDecayHalfAngle;
00287 static const G4ThreeVector origin;
00288
00289 G4int NSourceBin;
00290 G4double SBin[100];
00291 G4double SProfile[100];
00292 G4int NDecayBin;
00293 G4double DBin[100];
00294 G4double DProfile[100];
00295
00296 std::vector<G4String> LoadedNuclei;
00297 std::vector<G4String> ValidVolumes;
00298 bool isAllVolumesMode;
00299
00300 G4RadioactiveDecayRate theDecayRate;
00301 G4RadioactiveDecayRates theDecayRateVector;
00302 G4RadioactiveDecayRateVector theDecayRateTable;
00303 G4RadioactiveDecayRateTable theDecayRateTableVector;
00304
00305
00306 std::vector<G4RadioactivityTable*> theRadioactivityTables;
00307 G4int decayWindows[99];
00308 static const G4double levelTolerance;
00309
00310
00311 std::map<G4int, G4String> theUserRadioactiveDataFiles;
00312
00313
00314
00315 G4double fRemainderLifeTime;
00316 G4int verboseLevel;
00317
00318
00319
00320 G4ParticleChangeForRadDecay fParticleChangeForRadDecay;
00321
00322
00323 inline
00324 G4double AtRestGetPhysicalInteractionLength(const G4Track& track,
00325 G4ForceCondition* condition)
00326 {
00327 fRemainderLifeTime =
00328 G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(track, condition);
00329 return fRemainderLifeTime;
00330 }
00331
00332 inline
00333 G4VParticleChange* AtRestDoIt(const G4Track& theTrack,
00334 const G4Step& theStep)
00335 {return DecayIt(theTrack, theStep);}
00336
00337 inline
00338 G4VParticleChange* PostStepDoIt(const G4Track& theTrack,
00339 const G4Step& theStep)
00340 {return DecayIt(theTrack, theStep);}
00341 };
00342
00343 #endif
00344