#include <G4MiscLHEPBuilder_CHIPS.hh>
Public Member Functions | |
G4MiscLHEPBuilder_CHIPS () | |
virtual | ~G4MiscLHEPBuilder_CHIPS () |
void | Build () |
Private Attributes | |
G4AntiProtonInelasticProcess * | theAntiProtonInelastic |
G4LEAntiProtonInelastic * | theLEAntiProtonModel |
G4HEAntiProtonInelastic * | theHEAntiProtonModel |
G4QInelastic * | inelastic |
G4AntiNeutronInelasticProcess * | theAntiNeutronInelastic |
G4LEAntiNeutronInelastic * | theLEAntiNeutronModel |
G4HEAntiNeutronInelastic * | theHEAntiNeutronModel |
G4LambdaInelasticProcess * | theLambdaInelastic |
G4LELambdaInelastic * | theLELambdaModel |
G4HELambdaInelastic * | theHELambdaModel |
G4AntiLambdaInelasticProcess * | theAntiLambdaInelastic |
G4LEAntiLambdaInelastic * | theLEAntiLambdaModel |
G4HEAntiLambdaInelastic * | theHEAntiLambdaModel |
G4SigmaMinusInelasticProcess * | theSigmaMinusInelastic |
G4LESigmaMinusInelastic * | theLESigmaMinusModel |
G4HESigmaMinusInelastic * | theHESigmaMinusModel |
G4AntiSigmaMinusInelasticProcess * | theAntiSigmaMinusInelastic |
G4LEAntiSigmaMinusInelastic * | theLEAntiSigmaMinusModel |
G4HEAntiSigmaMinusInelastic * | theHEAntiSigmaMinusModel |
G4SigmaPlusInelasticProcess * | theSigmaPlusInelastic |
G4LESigmaPlusInelastic * | theLESigmaPlusModel |
G4HESigmaPlusInelastic * | theHESigmaPlusModel |
G4AntiSigmaPlusInelasticProcess * | theAntiSigmaPlusInelastic |
G4LEAntiSigmaPlusInelastic * | theLEAntiSigmaPlusModel |
G4HEAntiSigmaPlusInelastic * | theHEAntiSigmaPlusModel |
G4XiZeroInelasticProcess * | theXiZeroInelastic |
G4LEXiZeroInelastic * | theLEXiZeroModel |
G4HEXiZeroInelastic * | theHEXiZeroModel |
G4AntiXiZeroInelasticProcess * | theAntiXiZeroInelastic |
G4LEAntiXiZeroInelastic * | theLEAntiXiZeroModel |
G4HEAntiXiZeroInelastic * | theHEAntiXiZeroModel |
G4XiMinusInelasticProcess * | theXiMinusInelastic |
G4LEXiMinusInelastic * | theLEXiMinusModel |
G4HEXiMinusInelastic * | theHEXiMinusModel |
G4AntiXiMinusInelasticProcess * | theAntiXiMinusInelastic |
G4LEAntiXiMinusInelastic * | theLEAntiXiMinusModel |
G4HEAntiXiMinusInelastic * | theHEAntiXiMinusModel |
G4OmegaMinusInelasticProcess * | theOmegaMinusInelastic |
G4LEOmegaMinusInelastic * | theLEOmegaMinusModel |
G4HEOmegaMinusInelastic * | theHEOmegaMinusModel |
G4AntiOmegaMinusInelasticProcess * | theAntiOmegaMinusInelastic |
G4LEAntiOmegaMinusInelastic * | theLEAntiOmegaMinusModel |
G4HEAntiOmegaMinusInelastic * | theHEAntiOmegaMinusModel |
G4bool | wasActivated |
Definition at line 94 of file G4MiscLHEPBuilder_CHIPS.hh.
G4MiscLHEPBuilder_CHIPS::G4MiscLHEPBuilder_CHIPS | ( | ) |
G4MiscLHEPBuilder_CHIPS::~G4MiscLHEPBuilder_CHIPS | ( | ) | [virtual] |
Definition at line 47 of file G4MiscLHEPBuilder_CHIPS.cc.
References inelastic, and wasActivated.
00048 { 00049 if(wasActivated) delete inelastic; 00050 }
void G4MiscLHEPBuilder_CHIPS::Build | ( | ) |
Definition at line 52 of file G4MiscLHEPBuilder_CHIPS.cc.
References inelastic, theAntiLambdaInelastic, theAntiNeutronInelastic, theAntiOmegaMinusInelastic, theAntiSigmaMinusInelastic, theAntiSigmaPlusInelastic, theAntiXiMinusInelastic, theAntiXiZeroInelastic, theHEAntiLambdaModel, theHEAntiNeutronModel, theHEAntiOmegaMinusModel, theHEAntiSigmaMinusModel, theHEAntiSigmaPlusModel, theHEAntiXiMinusModel, theHEAntiXiZeroModel, theHELambdaModel, theHEOmegaMinusModel, theHESigmaMinusModel, theHESigmaPlusModel, theHEXiMinusModel, theHEXiZeroModel, theLambdaInelastic, theLEAntiLambdaModel, theLEAntiNeutronModel, theLEAntiOmegaMinusModel, theLEAntiSigmaMinusModel, theLEAntiSigmaPlusModel, theLEAntiXiMinusModel, theLEAntiXiZeroModel, theLELambdaModel, theLEOmegaMinusModel, theLESigmaMinusModel, theLESigmaPlusModel, theLEXiMinusModel, theLEXiZeroModel, theOmegaMinusInelastic, theSigmaMinusInelastic, theSigmaPlusInelastic, theXiMinusInelastic, theXiZeroInelastic, and wasActivated.
Referenced by HadronPhysicsQGSP_BERT_CHIPS::ConstructProcess().
00053 { 00054 G4ProcessManager * aProcMan = 0; 00055 wasActivated = true; 00056 00057 // anti-Proton 00058 //theAntiProtonInelastic = new G4AntiProtonInelasticProcess(); 00059 aProcMan = G4AntiProton::AntiProton()->GetProcessManager(); 00060 //theLEAntiProtonModel = new G4LEAntiProtonInelastic(); 00061 //theHEAntiProtonModel = new G4HEAntiProtonInelastic(); 00062 //theHEAntiProtonModel->SetMaxEnergy(100*TeV); 00063 //theAntiProtonInelastic->RegisterMe(theLEAntiProtonModel); 00064 //theAntiProtonInelastic->RegisterMe(theHEAntiProtonModel); 00065 //aProcMan->AddDiscreteProcess(theAntiProtonInelastic); 00066 //CHIPS 00067 inelastic = new G4QInelastic(); 00068 aProcMan->AddDiscreteProcess(inelastic); 00069 00070 // AntiNeutron 00071 theAntiNeutronInelastic = new G4AntiNeutronInelasticProcess(); 00072 aProcMan = G4AntiNeutron::AntiNeutron()->GetProcessManager(); 00073 theLEAntiNeutronModel = new G4LEAntiNeutronInelastic(); 00074 theHEAntiNeutronModel = new G4HEAntiNeutronInelastic(); 00075 theHEAntiNeutronModel->SetMaxEnergy(100*TeV); 00076 theAntiNeutronInelastic->RegisterMe(theLEAntiNeutronModel); 00077 theAntiNeutronInelastic->RegisterMe(theHEAntiNeutronModel); 00078 aProcMan->AddDiscreteProcess(theAntiNeutronInelastic); 00079 00080 // Lambda 00081 theLambdaInelastic = new G4LambdaInelasticProcess(); 00082 aProcMan = G4Lambda::Lambda()->GetProcessManager(); 00083 theLELambdaModel = new G4LELambdaInelastic(); 00084 theHELambdaModel = new G4HELambdaInelastic(); 00085 theHELambdaModel->SetMaxEnergy(100*TeV); 00086 theLambdaInelastic->RegisterMe(theLELambdaModel); 00087 theLambdaInelastic->RegisterMe(theHELambdaModel); 00088 aProcMan->AddDiscreteProcess(theLambdaInelastic); 00089 00090 // AntiLambda 00091 theAntiLambdaInelastic = new G4AntiLambdaInelasticProcess(); 00092 aProcMan = G4AntiLambda::AntiLambda()->GetProcessManager(); 00093 theLEAntiLambdaModel = new G4LEAntiLambdaInelastic(); 00094 theHEAntiLambdaModel = new G4HEAntiLambdaInelastic(); 00095 theHEAntiLambdaModel->SetMaxEnergy(100*TeV); 00096 theAntiLambdaInelastic->RegisterMe(theLEAntiLambdaModel); 00097 theAntiLambdaInelastic->RegisterMe(theHEAntiLambdaModel); 00098 aProcMan->AddDiscreteProcess(theAntiLambdaInelastic); 00099 00100 // SigmaMinus 00101 theSigmaMinusInelastic = new G4SigmaMinusInelasticProcess(); 00102 aProcMan = G4SigmaMinus::SigmaMinus()->GetProcessManager(); 00103 theLESigmaMinusModel = new G4LESigmaMinusInelastic(); 00104 theHESigmaMinusModel = new G4HESigmaMinusInelastic(); 00105 theHESigmaMinusModel->SetMaxEnergy(100*TeV); 00106 theSigmaMinusInelastic->RegisterMe(theLESigmaMinusModel); 00107 theSigmaMinusInelastic->RegisterMe(theHESigmaMinusModel); 00108 aProcMan->AddDiscreteProcess(theSigmaMinusInelastic); 00109 00110 // anti-SigmaMinus 00111 theAntiSigmaMinusInelastic = new G4AntiSigmaMinusInelasticProcess(); 00112 aProcMan = G4AntiSigmaMinus::AntiSigmaMinus()->GetProcessManager(); 00113 theLEAntiSigmaMinusModel = new G4LEAntiSigmaMinusInelastic(); 00114 theHEAntiSigmaMinusModel = new G4HEAntiSigmaMinusInelastic(); 00115 theHEAntiSigmaMinusModel->SetMaxEnergy(100*TeV); 00116 theAntiSigmaMinusInelastic->RegisterMe(theLEAntiSigmaMinusModel); 00117 theAntiSigmaMinusInelastic->RegisterMe(theHEAntiSigmaMinusModel); 00118 aProcMan->AddDiscreteProcess(theAntiSigmaMinusInelastic); 00119 00120 // SigmaPlus 00121 theSigmaPlusInelastic = new G4SigmaPlusInelasticProcess(); 00122 aProcMan = G4SigmaPlus::SigmaPlus()->GetProcessManager(); 00123 theLESigmaPlusModel = new G4LESigmaPlusInelastic(); 00124 theHESigmaPlusModel = new G4HESigmaPlusInelastic(); 00125 theHESigmaPlusModel->SetMaxEnergy(100*TeV); 00126 theSigmaPlusInelastic->RegisterMe(theLESigmaPlusModel); 00127 theSigmaPlusInelastic->RegisterMe(theHESigmaPlusModel); 00128 aProcMan->AddDiscreteProcess(theSigmaPlusInelastic); 00129 00130 // anti-SigmaPlus 00131 theAntiSigmaPlusInelastic = new G4AntiSigmaPlusInelasticProcess(); 00132 aProcMan = G4AntiSigmaPlus::AntiSigmaPlus()->GetProcessManager(); 00133 theLEAntiSigmaPlusModel = new G4LEAntiSigmaPlusInelastic(); 00134 theHEAntiSigmaPlusModel = new G4HEAntiSigmaPlusInelastic(); 00135 theHEAntiSigmaPlusModel->SetMaxEnergy(100*TeV); 00136 theAntiSigmaPlusInelastic->RegisterMe(theLEAntiSigmaPlusModel); 00137 theAntiSigmaPlusInelastic->RegisterMe(theHEAntiSigmaPlusModel); 00138 aProcMan->AddDiscreteProcess(theAntiSigmaPlusInelastic); 00139 00140 // XiMinus 00141 theXiMinusInelastic = new G4XiMinusInelasticProcess(); 00142 aProcMan = G4XiMinus::XiMinus()->GetProcessManager(); 00143 theLEXiMinusModel = new G4LEXiMinusInelastic(); 00144 theHEXiMinusModel = new G4HEXiMinusInelastic(); 00145 theHEXiMinusModel->SetMaxEnergy(100*TeV); 00146 theXiMinusInelastic->RegisterMe(theLEXiMinusModel); 00147 theXiMinusInelastic->RegisterMe(theHEXiMinusModel); 00148 aProcMan->AddDiscreteProcess(theXiMinusInelastic); 00149 00150 // anti-XiMinus 00151 theAntiXiMinusInelastic = new G4AntiXiMinusInelasticProcess(); 00152 aProcMan = G4AntiXiMinus::AntiXiMinus()->GetProcessManager(); 00153 theLEAntiXiMinusModel = new G4LEAntiXiMinusInelastic(); 00154 theHEAntiXiMinusModel = new G4HEAntiXiMinusInelastic(); 00155 theHEAntiXiMinusModel->SetMaxEnergy(100*TeV); 00156 theAntiXiMinusInelastic->RegisterMe(theLEAntiXiMinusModel); 00157 theAntiXiMinusInelastic->RegisterMe(theHEAntiXiMinusModel); 00158 aProcMan->AddDiscreteProcess(theAntiXiMinusInelastic); 00159 00160 // XiZero 00161 theXiZeroInelastic = new G4XiZeroInelasticProcess(); 00162 aProcMan = G4XiZero::XiZero()->GetProcessManager(); 00163 theLEXiZeroModel = new G4LEXiZeroInelastic(); 00164 theHEXiZeroModel = new G4HEXiZeroInelastic(); 00165 theHEXiZeroModel->SetMaxEnergy(100*TeV); 00166 theXiZeroInelastic->RegisterMe(theLEXiZeroModel); 00167 theXiZeroInelastic->RegisterMe(theHEXiZeroModel); 00168 aProcMan->AddDiscreteProcess(theXiZeroInelastic); 00169 00170 // anti-XiZero 00171 theAntiXiZeroInelastic = new G4AntiXiZeroInelasticProcess(); 00172 aProcMan = G4AntiXiZero::AntiXiZero()->GetProcessManager(); 00173 theLEAntiXiZeroModel = new G4LEAntiXiZeroInelastic(); 00174 theHEAntiXiZeroModel = new G4HEAntiXiZeroInelastic(); 00175 theHEAntiXiZeroModel->SetMaxEnergy(100*TeV); 00176 theAntiXiZeroInelastic->RegisterMe(theLEAntiXiZeroModel); 00177 theAntiXiZeroInelastic->RegisterMe(theHEAntiXiZeroModel); 00178 aProcMan->AddDiscreteProcess(theAntiXiZeroInelastic); 00179 00180 // OmegaMinus 00181 theOmegaMinusInelastic = new G4OmegaMinusInelasticProcess(); 00182 aProcMan = G4OmegaMinus::OmegaMinus()->GetProcessManager(); 00183 theLEOmegaMinusModel = new G4LEOmegaMinusInelastic(); 00184 theHEOmegaMinusModel = new G4HEOmegaMinusInelastic(); 00185 theHEOmegaMinusModel->SetMaxEnergy(100*TeV); 00186 theOmegaMinusInelastic->RegisterMe(theLEOmegaMinusModel); 00187 theOmegaMinusInelastic->RegisterMe(theHEOmegaMinusModel); 00188 aProcMan->AddDiscreteProcess(theOmegaMinusInelastic); 00189 00190 // anti-OmegaMinus 00191 theAntiOmegaMinusInelastic = new G4AntiOmegaMinusInelasticProcess(); 00192 aProcMan = G4AntiOmegaMinus::AntiOmegaMinus()->GetProcessManager(); 00193 theLEAntiOmegaMinusModel = new G4LEAntiOmegaMinusInelastic(); 00194 theHEAntiOmegaMinusModel = new G4HEAntiOmegaMinusInelastic(); 00195 theHEAntiOmegaMinusModel->SetMaxEnergy(100*TeV); 00196 theAntiOmegaMinusInelastic->RegisterMe(theLEAntiOmegaMinusModel); 00197 theAntiOmegaMinusInelastic->RegisterMe(theHEAntiOmegaMinusModel); 00198 aProcMan->AddDiscreteProcess(theAntiOmegaMinusInelastic); 00199 }
G4QInelastic* G4MiscLHEPBuilder_CHIPS::inelastic [private] |
Definition at line 109 of file G4MiscLHEPBuilder_CHIPS.hh.
Referenced by Build(), and ~G4MiscLHEPBuilder_CHIPS().
G4AntiLambdaInelasticProcess* G4MiscLHEPBuilder_CHIPS::theAntiLambdaInelastic [private] |
G4AntiNeutronInelasticProcess* G4MiscLHEPBuilder_CHIPS::theAntiNeutronInelastic [private] |
G4AntiOmegaMinusInelasticProcess* G4MiscLHEPBuilder_CHIPS::theAntiOmegaMinusInelastic [private] |
G4AntiProtonInelasticProcess* G4MiscLHEPBuilder_CHIPS::theAntiProtonInelastic [private] |
Definition at line 106 of file G4MiscLHEPBuilder_CHIPS.hh.
G4AntiSigmaMinusInelasticProcess* G4MiscLHEPBuilder_CHIPS::theAntiSigmaMinusInelastic [private] |
G4AntiSigmaPlusInelasticProcess* G4MiscLHEPBuilder_CHIPS::theAntiSigmaPlusInelastic [private] |
G4AntiXiMinusInelasticProcess* G4MiscLHEPBuilder_CHIPS::theAntiXiMinusInelastic [private] |
G4AntiXiZeroInelasticProcess* G4MiscLHEPBuilder_CHIPS::theAntiXiZeroInelastic [private] |
G4HEAntiLambdaInelastic* G4MiscLHEPBuilder_CHIPS::theHEAntiLambdaModel [private] |
G4HEAntiNeutronInelastic* G4MiscLHEPBuilder_CHIPS::theHEAntiNeutronModel [private] |
G4HEAntiOmegaMinusInelastic* G4MiscLHEPBuilder_CHIPS::theHEAntiOmegaMinusModel [private] |
G4HEAntiProtonInelastic* G4MiscLHEPBuilder_CHIPS::theHEAntiProtonModel [private] |
Definition at line 108 of file G4MiscLHEPBuilder_CHIPS.hh.
G4HEAntiSigmaMinusInelastic* G4MiscLHEPBuilder_CHIPS::theHEAntiSigmaMinusModel [private] |
G4HEAntiSigmaPlusInelastic* G4MiscLHEPBuilder_CHIPS::theHEAntiSigmaPlusModel [private] |
G4HEAntiXiMinusInelastic* G4MiscLHEPBuilder_CHIPS::theHEAntiXiMinusModel [private] |
G4HEAntiXiZeroInelastic* G4MiscLHEPBuilder_CHIPS::theHEAntiXiZeroModel [private] |
G4HELambdaInelastic* G4MiscLHEPBuilder_CHIPS::theHELambdaModel [private] |
G4HEOmegaMinusInelastic* G4MiscLHEPBuilder_CHIPS::theHEOmegaMinusModel [private] |
G4HESigmaMinusInelastic* G4MiscLHEPBuilder_CHIPS::theHESigmaMinusModel [private] |
G4HESigmaPlusInelastic* G4MiscLHEPBuilder_CHIPS::theHESigmaPlusModel [private] |
G4HEXiMinusInelastic* G4MiscLHEPBuilder_CHIPS::theHEXiMinusModel [private] |
G4HEXiZeroInelastic* G4MiscLHEPBuilder_CHIPS::theHEXiZeroModel [private] |
G4LambdaInelasticProcess* G4MiscLHEPBuilder_CHIPS::theLambdaInelastic [private] |
G4LEAntiLambdaInelastic* G4MiscLHEPBuilder_CHIPS::theLEAntiLambdaModel [private] |
G4LEAntiNeutronInelastic* G4MiscLHEPBuilder_CHIPS::theLEAntiNeutronModel [private] |
G4LEAntiOmegaMinusInelastic* G4MiscLHEPBuilder_CHIPS::theLEAntiOmegaMinusModel [private] |
G4LEAntiProtonInelastic* G4MiscLHEPBuilder_CHIPS::theLEAntiProtonModel [private] |
Definition at line 107 of file G4MiscLHEPBuilder_CHIPS.hh.
G4LEAntiSigmaMinusInelastic* G4MiscLHEPBuilder_CHIPS::theLEAntiSigmaMinusModel [private] |
G4LEAntiSigmaPlusInelastic* G4MiscLHEPBuilder_CHIPS::theLEAntiSigmaPlusModel [private] |
G4LEAntiXiMinusInelastic* G4MiscLHEPBuilder_CHIPS::theLEAntiXiMinusModel [private] |
G4LEAntiXiZeroInelastic* G4MiscLHEPBuilder_CHIPS::theLEAntiXiZeroModel [private] |
G4LELambdaInelastic* G4MiscLHEPBuilder_CHIPS::theLELambdaModel [private] |
G4LEOmegaMinusInelastic* G4MiscLHEPBuilder_CHIPS::theLEOmegaMinusModel [private] |
G4LESigmaMinusInelastic* G4MiscLHEPBuilder_CHIPS::theLESigmaMinusModel [private] |
G4LESigmaPlusInelastic* G4MiscLHEPBuilder_CHIPS::theLESigmaPlusModel [private] |
G4LEXiMinusInelastic* G4MiscLHEPBuilder_CHIPS::theLEXiMinusModel [private] |
G4LEXiZeroInelastic* G4MiscLHEPBuilder_CHIPS::theLEXiZeroModel [private] |
G4OmegaMinusInelasticProcess* G4MiscLHEPBuilder_CHIPS::theOmegaMinusInelastic [private] |
G4SigmaMinusInelasticProcess* G4MiscLHEPBuilder_CHIPS::theSigmaMinusInelastic [private] |
G4SigmaPlusInelasticProcess* G4MiscLHEPBuilder_CHIPS::theSigmaPlusInelastic [private] |
G4XiMinusInelasticProcess* G4MiscLHEPBuilder_CHIPS::theXiMinusInelastic [private] |
G4XiZeroInelasticProcess* G4MiscLHEPBuilder_CHIPS::theXiZeroInelastic [private] |
G4bool G4MiscLHEPBuilder_CHIPS::wasActivated [private] |
Definition at line 176 of file G4MiscLHEPBuilder_CHIPS.hh.
Referenced by Build(), and ~G4MiscLHEPBuilder_CHIPS().