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
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048 #ifndef G4HadronicProcess_h
00049 #define G4HadronicProcess_h 1
00050
00051 #include "globals.hh"
00052 #include "G4VDiscreteProcess.hh"
00053 #include "G4EnergyRangeManager.hh"
00054 #include "G4Nucleus.hh"
00055 #include "G4ReactionProduct.hh"
00056 #include <vector>
00057 #include "G4VCrossSectionDataSet.hh"
00058 #include "G4VLeadingParticleBiasing.hh"
00059
00060 #include "G4CrossSectionDataStore.hh"
00061 #include "G4HadronicProcessType.hh"
00062
00063 class G4Track;
00064 class G4Step;
00065 class G4Element;
00066 class G4ParticleChange;
00067
00068
00069 class G4HadronicProcess : public G4VDiscreteProcess
00070 {
00071 public:
00072 G4HadronicProcess(const G4String& processName="Hadronic",
00073 G4ProcessType procType=fHadronic);
00074
00075
00076 G4HadronicProcess(const G4String& processName,
00077 G4HadronicProcessType subType);
00078
00079 virtual ~G4HadronicProcess();
00080
00081
00082 void RegisterMe(G4HadronicInteraction* a);
00083
00084
00085 inline
00086 G4double GetElementCrossSection(const G4DynamicParticle * part,
00087 const G4Element * elm,
00088 const G4Material* mat = 0)
00089 {
00090 G4double x = theCrossSectionDataStore->GetCrossSection(part, elm, mat);
00091 if(x < 0.0) { x = 0.0; }
00092 return x;
00093 }
00094
00095
00096 inline
00097 G4double GetMicroscopicCrossSection(const G4DynamicParticle * part,
00098 const G4Element * elm,
00099 const G4Material* mat = 0)
00100 { return GetElementCrossSection(part, elm, mat); }
00101
00102
00103 virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
00104 const G4Step& aStep);
00105
00106
00107 virtual void PreparePhysicsTable(const G4ParticleDefinition&);
00108
00109
00110 virtual void BuildPhysicsTable(const G4ParticleDefinition&);
00111
00112
00113 inline void DumpPhysicsTable(const G4ParticleDefinition& p)
00114 { theCrossSectionDataStore->DumpPhysicsTable(p); }
00115
00116
00117 inline void AddDataSet(G4VCrossSectionDataSet * aDataSet)
00118 { theCrossSectionDataStore->AddDataSet(aDataSet);}
00119
00120
00121 inline G4EnergyRangeManager *GetManagerPointer()
00122 { return &theEnergyRangeManager; }
00123
00124
00125 G4double GetMeanFreePath(const G4Track &aTrack, G4double,
00126 G4ForceCondition *);
00127
00128
00129 inline const G4Nucleus* GetTargetNucleus() const
00130 { return &targetNucleus; }
00131
00132
00133 inline const G4Isotope* GetTargetIsotope()
00134 { return targetNucleus.GetIsotope(); }
00135
00136 virtual void ProcessDescription(std::ostream& outFile) const;
00137
00138 protected:
00139
00140
00141
00142 inline G4HadronicInteraction* ChooseHadronicInteraction(
00143 G4double kineticEnergy, G4Material* aMaterial, G4Element* anElement)
00144 { return theEnergyRangeManager.GetHadronicInteraction(kineticEnergy,
00145 aMaterial,anElement);
00146 }
00147
00148
00149 inline G4Nucleus* GetTargetNucleusPointer()
00150 { return &targetNucleus; }
00151
00152 public:
00153
00154 void BiasCrossSectionByFactor(G4double aScale);
00155
00156
00157 inline void SetEpReportLevel(G4int level)
00158 { epReportLevel = level; }
00159
00160 inline void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
00161 { epCheckLevels.first = relativeLevel;
00162 epCheckLevels.second = absoluteLevel;
00163 levelsSetByProcess = true;
00164 }
00165
00166 inline std::pair<G4double, G4double> GetEnergyMomentumCheckLevels() const
00167 { return epCheckLevels; }
00168
00169
00170 inline G4CrossSectionDataStore* GetCrossSectionDataStore()
00171 {return theCrossSectionDataStore;}
00172
00173 inline void MultiplyCrossSectionBy(G4double factor)
00174 { aScaleFactor = factor; }
00175
00176 protected:
00177
00178 void DumpState(const G4Track&, const G4String&, G4ExceptionDescription&);
00179
00180
00181 inline const G4EnergyRangeManager &GetEnergyRangeManager() const
00182 { return theEnergyRangeManager; }
00183
00184
00185 inline void SetEnergyRangeManager( const G4EnergyRangeManager &value )
00186 { theEnergyRangeManager = value; }
00187
00188
00189 inline G4HadronicInteraction* GetHadronicInteraction() const
00190 { return theInteraction; }
00191
00192
00193 inline G4double GetLastCrossSection()
00194 { return theLastCrossSection; }
00195
00196
00197 void FillResult(G4HadFinalState* aR, const G4Track& aT);
00198
00199
00200 G4HadFinalState* CheckResult(const G4HadProjectile& thePro,
00201 const G4Nucleus& targetNucleus,
00202 G4HadFinalState* result) const;
00203
00204
00205 void CheckEnergyMomentumConservation(const G4Track&, const G4Nucleus&);
00206
00207 private:
00208 G4double XBiasSurvivalProbability();
00209 G4double XBiasSecondaryWeight();
00210
00211
00212 G4HadronicProcess& operator=(const G4HadronicProcess& right);
00213 G4HadronicProcess(const G4HadronicProcess&);
00214
00215
00216 void GetEnergyMomentumCheckEnvvars();
00217
00218 protected:
00219
00220 G4HadProjectile thePro;
00221
00222 G4ParticleChange* theTotalResult;
00223
00224 G4int epReportLevel;
00225
00226 private:
00227
00228 G4EnergyRangeManager theEnergyRangeManager;
00229
00230 G4HadronicInteraction* theInteraction;
00231
00232 G4CrossSectionDataStore* theCrossSectionDataStore;
00233
00234 G4Nucleus targetNucleus;
00235
00236 bool G4HadronicProcess_debug_flag;
00237
00238
00239 std::pair<G4double, G4double> epCheckLevels;
00240 G4bool levelsSetByProcess;
00241
00242 std::vector<G4VLeadingParticleBiasing *> theBias;
00243
00244 G4double theInitialNumberOfInteractionLength;
00245
00246 G4double aScaleFactor;
00247 G4bool xBiasOn;
00248 G4double theLastCrossSection;
00249 };
00250
00251 #endif
00252