#include <G4MiscBuilder.hh>
Public Member Functions | |
G4MiscBuilder () | |
virtual | ~G4MiscBuilder () |
void | Build () |
Definition at line 96 of file G4MiscBuilder.hh.
G4MiscBuilder::G4MiscBuilder | ( | ) |
Definition at line 47 of file G4MiscBuilder.cc.
00047 : 00048 theAntiProtonInelastic(0), theLEAntiProtonModel(0), 00049 theHEAntiProtonModel(0), 00050 theAntiNeutronInelastic(0), theLEAntiNeutronModel(0), 00051 theHEAntiNeutronModel(0), 00052 theLambdaInelastic(0), theLELambdaModel(0), theHELambdaModel(0), 00053 theAntiLambdaInelastic(0), theLEAntiLambdaModel(0), theHEAntiLambdaModel(0), 00054 theSigmaMinusInelastic(0), theLESigmaMinusModel(0), theHESigmaMinusModel(0), 00055 theAntiSigmaMinusInelastic(0), theLEAntiSigmaMinusModel(0), theHEAntiSigmaMinusModel(0), 00056 theSigmaPlusInelastic(0), theLESigmaPlusModel(0), theHESigmaPlusModel(0), 00057 theAntiSigmaPlusInelastic(0), theLEAntiSigmaPlusModel(0), theHEAntiSigmaPlusModel(0), 00058 theXiZeroInelastic(0), theLEXiZeroModel(0), theHEXiZeroModel(0), 00059 theAntiXiZeroInelastic(0), theLEAntiXiZeroModel(0), theHEAntiXiZeroModel(0), 00060 theXiMinusInelastic(0), theLEXiMinusModel(0), theHEXiMinusModel(0), 00061 theAntiXiMinusInelastic(0), theLEAntiXiMinusModel(0), theHEAntiXiMinusModel(0), 00062 theOmegaMinusInelastic(0), theLEOmegaMinusModel(0), theHEOmegaMinusModel(0), 00063 theAntiOmegaMinusInelastic(0), theLEAntiOmegaMinusModel(0), theHEAntiOmegaMinusModel(0), 00064 wasActivated(false) 00065 00066 {}
G4MiscBuilder::~G4MiscBuilder | ( | ) | [virtual] |
void G4MiscBuilder::Build | ( | ) |
Definition at line 73 of file G4MiscBuilder.cc.
References G4HadronicProcess::AddDataSet(), G4ProcessManager::AddDiscreteProcess(), G4AntiLambda::AntiLambda(), G4AntiNeutron::AntiNeutron(), G4AntiOmegaMinus::AntiOmegaMinus(), G4AntiProton::AntiProton(), G4AntiSigmaMinus::AntiSigmaMinus(), G4AntiSigmaPlus::AntiSigmaPlus(), G4AntiXiMinus::AntiXiMinus(), G4AntiXiZero::AntiXiZero(), G4ChipsHyperonInelasticXS::Default_Name(), G4CrossSectionDataSetRegistry::GetCrossSectionDataSet(), G4ParticleDefinition::GetProcessManager(), G4CrossSectionDataSetRegistry::Instance(), G4Lambda::Lambda(), G4OmegaMinus::OmegaMinus(), G4HadronicProcess::RegisterMe(), G4HadronicInteraction::SetMaxEnergy(), G4SigmaMinus::SigmaMinus(), G4SigmaPlus::SigmaPlus(), G4XiMinus::XiMinus(), and G4XiZero::XiZero().
Referenced by HadronPhysicsQGSP_INCLXX::ConstructProcess(), HadronPhysicsQGSP_BIC_HP::ConstructProcess(), HadronPhysicsQGSP_BIC::ConstructProcess(), HadronPhysicsQGSP_BERT_TRV::ConstructProcess(), HadronPhysicsQGSP_BERT_NOLEP::ConstructProcess(), HadronPhysicsQGSP_BERT_HP::ConstructProcess(), HadronPhysicsQGSP_BERT::ConstructProcess(), HadronPhysicsQGSP::ConstructProcess(), and HadronPhysicsQGS_BIC::ConstructProcess().
00074 { 00075 G4ProcessManager * aProcMan = 0; 00076 wasActivated = true; 00077 00078 theAntiNucleonData = new G4CrossSectionInelastic(new G4ComponentAntiNuclNuclearXS()); 00079 theChipsInelastic = G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4ChipsHyperonInelasticXS::Default_Name()); 00080 00081 // anti-Proton 00082 theAntiProtonInelastic = new G4AntiProtonInelasticProcess(); 00083 aProcMan = G4AntiProton::AntiProton()->GetProcessManager(); 00084 theLEAntiProtonModel = new G4LEAntiProtonInelastic(); 00085 theHEAntiProtonModel = new G4HEAntiProtonInelastic(); 00086 theHEAntiProtonModel->SetMaxEnergy(100*TeV); 00087 theAntiProtonInelastic->RegisterMe(theLEAntiProtonModel); 00088 theAntiProtonInelastic->RegisterMe(theHEAntiProtonModel); 00089 aProcMan->AddDiscreteProcess(theAntiProtonInelastic); 00090 theAntiProtonInelastic->AddDataSet(theAntiNucleonData); 00091 00092 // AntiNeutron 00093 theAntiNeutronInelastic = new G4AntiNeutronInelasticProcess(); 00094 aProcMan = G4AntiNeutron::AntiNeutron()->GetProcessManager(); 00095 theLEAntiNeutronModel = new G4LEAntiNeutronInelastic(); 00096 theHEAntiNeutronModel = new G4HEAntiNeutronInelastic(); 00097 theHEAntiNeutronModel->SetMaxEnergy(100*TeV); 00098 theAntiNeutronInelastic->RegisterMe(theLEAntiNeutronModel); 00099 theAntiNeutronInelastic->RegisterMe(theHEAntiNeutronModel); 00100 aProcMan->AddDiscreteProcess(theAntiNeutronInelastic); 00101 theAntiNeutronInelastic->AddDataSet(theAntiNucleonData); 00102 00103 // Lambda 00104 theLambdaInelastic = new G4LambdaInelasticProcess(); 00105 aProcMan = G4Lambda::Lambda()->GetProcessManager(); 00106 theLELambdaModel = new G4LELambdaInelastic(); 00107 theHELambdaModel = new G4HELambdaInelastic(); 00108 theHELambdaModel->SetMaxEnergy(100*TeV); 00109 theLambdaInelastic->RegisterMe(theLELambdaModel); 00110 theLambdaInelastic->RegisterMe(theHELambdaModel); 00111 aProcMan->AddDiscreteProcess(theLambdaInelastic); 00112 theLambdaInelastic->AddDataSet(theChipsInelastic); 00113 00114 // AntiLambda 00115 theAntiLambdaInelastic = new G4AntiLambdaInelasticProcess(); 00116 aProcMan = G4AntiLambda::AntiLambda()->GetProcessManager(); 00117 theLEAntiLambdaModel = new G4LEAntiLambdaInelastic(); 00118 theHEAntiLambdaModel = new G4HEAntiLambdaInelastic(); 00119 theHEAntiLambdaModel->SetMaxEnergy(100*TeV); 00120 theAntiLambdaInelastic->RegisterMe(theLEAntiLambdaModel); 00121 theAntiLambdaInelastic->RegisterMe(theHEAntiLambdaModel); 00122 aProcMan->AddDiscreteProcess(theAntiLambdaInelastic); 00123 theAntiLambdaInelastic->AddDataSet(theChipsInelastic); 00124 00125 // SigmaMinus 00126 theSigmaMinusInelastic = new G4SigmaMinusInelasticProcess(); 00127 aProcMan = G4SigmaMinus::SigmaMinus()->GetProcessManager(); 00128 theLESigmaMinusModel = new G4LESigmaMinusInelastic(); 00129 theHESigmaMinusModel = new G4HESigmaMinusInelastic(); 00130 theHESigmaMinusModel->SetMaxEnergy(100*TeV); 00131 theSigmaMinusInelastic->RegisterMe(theLESigmaMinusModel); 00132 theSigmaMinusInelastic->RegisterMe(theHESigmaMinusModel); 00133 aProcMan->AddDiscreteProcess(theSigmaMinusInelastic); 00134 theSigmaMinusInelastic->AddDataSet(theChipsInelastic); 00135 00136 // anti-SigmaMinus 00137 theAntiSigmaMinusInelastic = new G4AntiSigmaMinusInelasticProcess(); 00138 aProcMan = G4AntiSigmaMinus::AntiSigmaMinus()->GetProcessManager(); 00139 theLEAntiSigmaMinusModel = new G4LEAntiSigmaMinusInelastic(); 00140 theHEAntiSigmaMinusModel = new G4HEAntiSigmaMinusInelastic(); 00141 theHEAntiSigmaMinusModel->SetMaxEnergy(100*TeV); 00142 theAntiSigmaMinusInelastic->RegisterMe(theLEAntiSigmaMinusModel); 00143 theAntiSigmaMinusInelastic->RegisterMe(theHEAntiSigmaMinusModel); 00144 aProcMan->AddDiscreteProcess(theAntiSigmaMinusInelastic); 00145 theAntiSigmaMinusInelastic->AddDataSet(theChipsInelastic); 00146 00147 // SigmaPlus 00148 theSigmaPlusInelastic = new G4SigmaPlusInelasticProcess(); 00149 aProcMan = G4SigmaPlus::SigmaPlus()->GetProcessManager(); 00150 theLESigmaPlusModel = new G4LESigmaPlusInelastic(); 00151 theHESigmaPlusModel = new G4HESigmaPlusInelastic(); 00152 theHESigmaPlusModel->SetMaxEnergy(100*TeV); 00153 theSigmaPlusInelastic->RegisterMe(theLESigmaPlusModel); 00154 theSigmaPlusInelastic->RegisterMe(theHESigmaPlusModel); 00155 aProcMan->AddDiscreteProcess(theSigmaPlusInelastic); 00156 theSigmaPlusInelastic->AddDataSet(theChipsInelastic); 00157 00158 // anti-SigmaPlus 00159 theAntiSigmaPlusInelastic = new G4AntiSigmaPlusInelasticProcess(); 00160 aProcMan = G4AntiSigmaPlus::AntiSigmaPlus()->GetProcessManager(); 00161 theLEAntiSigmaPlusModel = new G4LEAntiSigmaPlusInelastic(); 00162 theHEAntiSigmaPlusModel = new G4HEAntiSigmaPlusInelastic(); 00163 theHEAntiSigmaPlusModel->SetMaxEnergy(100*TeV); 00164 theAntiSigmaPlusInelastic->RegisterMe(theLEAntiSigmaPlusModel); 00165 theAntiSigmaPlusInelastic->RegisterMe(theHEAntiSigmaPlusModel); 00166 aProcMan->AddDiscreteProcess(theAntiSigmaPlusInelastic); 00167 theAntiSigmaPlusInelastic->AddDataSet(theChipsInelastic); 00168 00169 // XiMinus 00170 theXiMinusInelastic = new G4XiMinusInelasticProcess(); 00171 aProcMan = G4XiMinus::XiMinus()->GetProcessManager(); 00172 theLEXiMinusModel = new G4LEXiMinusInelastic(); 00173 theHEXiMinusModel = new G4HEXiMinusInelastic(); 00174 theHEXiMinusModel->SetMaxEnergy(100*TeV); 00175 theXiMinusInelastic->RegisterMe(theLEXiMinusModel); 00176 theXiMinusInelastic->RegisterMe(theHEXiMinusModel); 00177 aProcMan->AddDiscreteProcess(theXiMinusInelastic); 00178 theXiMinusInelastic->AddDataSet(theChipsInelastic); 00179 00180 // anti-XiMinus 00181 theAntiXiMinusInelastic = new G4AntiXiMinusInelasticProcess(); 00182 aProcMan = G4AntiXiMinus::AntiXiMinus()->GetProcessManager(); 00183 theLEAntiXiMinusModel = new G4LEAntiXiMinusInelastic(); 00184 theHEAntiXiMinusModel = new G4HEAntiXiMinusInelastic(); 00185 theHEAntiXiMinusModel->SetMaxEnergy(100*TeV); 00186 theAntiXiMinusInelastic->RegisterMe(theLEAntiXiMinusModel); 00187 theAntiXiMinusInelastic->RegisterMe(theHEAntiXiMinusModel); 00188 aProcMan->AddDiscreteProcess(theAntiXiMinusInelastic); 00189 theAntiXiMinusInelastic->AddDataSet(theChipsInelastic); 00190 00191 // XiZero 00192 theXiZeroInelastic = new G4XiZeroInelasticProcess(); 00193 aProcMan = G4XiZero::XiZero()->GetProcessManager(); 00194 theLEXiZeroModel = new G4LEXiZeroInelastic(); 00195 theHEXiZeroModel = new G4HEXiZeroInelastic(); 00196 theHEXiZeroModel->SetMaxEnergy(100*TeV); 00197 theXiZeroInelastic->RegisterMe(theLEXiZeroModel); 00198 theXiZeroInelastic->RegisterMe(theHEXiZeroModel); 00199 aProcMan->AddDiscreteProcess(theXiZeroInelastic); 00200 theXiZeroInelastic->AddDataSet(theChipsInelastic); 00201 00202 // anti-XiZero 00203 theAntiXiZeroInelastic = new G4AntiXiZeroInelasticProcess(); 00204 aProcMan = G4AntiXiZero::AntiXiZero()->GetProcessManager(); 00205 theLEAntiXiZeroModel = new G4LEAntiXiZeroInelastic(); 00206 theHEAntiXiZeroModel = new G4HEAntiXiZeroInelastic(); 00207 theHEAntiXiZeroModel->SetMaxEnergy(100*TeV); 00208 theAntiXiZeroInelastic->RegisterMe(theLEAntiXiZeroModel); 00209 theAntiXiZeroInelastic->RegisterMe(theHEAntiXiZeroModel); 00210 aProcMan->AddDiscreteProcess(theAntiXiZeroInelastic); 00211 theAntiXiZeroInelastic->AddDataSet(theChipsInelastic); 00212 00213 // OmegaMinus 00214 theOmegaMinusInelastic = new G4OmegaMinusInelasticProcess(); 00215 aProcMan = G4OmegaMinus::OmegaMinus()->GetProcessManager(); 00216 theLEOmegaMinusModel = new G4LEOmegaMinusInelastic(); 00217 theHEOmegaMinusModel = new G4HEOmegaMinusInelastic(); 00218 theHEOmegaMinusModel->SetMaxEnergy(100*TeV); 00219 theOmegaMinusInelastic->RegisterMe(theLEOmegaMinusModel); 00220 theOmegaMinusInelastic->RegisterMe(theHEOmegaMinusModel); 00221 aProcMan->AddDiscreteProcess(theOmegaMinusInelastic); 00222 theOmegaMinusInelastic->AddDataSet(theChipsInelastic); 00223 00224 // anti-OmegaMinus 00225 theAntiOmegaMinusInelastic = new G4AntiOmegaMinusInelasticProcess(); 00226 aProcMan = G4AntiOmegaMinus::AntiOmegaMinus()->GetProcessManager(); 00227 theLEAntiOmegaMinusModel = new G4LEAntiOmegaMinusInelastic(); 00228 theHEAntiOmegaMinusModel = new G4HEAntiOmegaMinusInelastic(); 00229 theHEAntiOmegaMinusModel->SetMaxEnergy(100*TeV); 00230 theAntiOmegaMinusInelastic->RegisterMe(theLEAntiOmegaMinusModel); 00231 theAntiOmegaMinusInelastic->RegisterMe(theHEAntiOmegaMinusModel); 00232 aProcMan->AddDiscreteProcess(theAntiOmegaMinusInelastic); 00233 theAntiOmegaMinusInelastic->AddDataSet(theChipsInelastic); 00234 }