00001 #include"G4ParticleTable.hh"
00002 #include "Randomize.hh"
00003
00004 #include<iostream>
00005 #include<fstream>
00006 #include <string>
00007
00008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00009 #include "FWCore/ParameterSet/interface/FileInPath.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011
00012 #include"SimG4Core/CustomPhysics/interface/G4ProcessHelper.hh"
00013 #include "SimG4Core/CustomPhysics/interface/CustomPDGParser.h"
00014
00015 G4ProcessHelper::G4ProcessHelper(const edm::ParameterSet & p){
00016
00017 particleTable = G4ParticleTable::GetParticleTable();
00018
00019 theProton = particleTable->FindParticle("proton");
00020 theNeutron = particleTable->FindParticle("neutron");
00021
00022 G4String line;
00023
00024 edm::FileInPath fp = p.getParameter<edm::FileInPath>("processesDef");
00025 std::string processDefFilePath = fp.fullPath();
00026 std::ifstream process_stream (processDefFilePath.c_str());
00027
00028 resonant = p.getParameter<bool>("resonant");
00029 ek_0 = p.getParameter<double>("resonanceEnergy")*GeV;
00030 gamma = p.getParameter<double>("gamma")*GeV;
00031 amplitude = p.getParameter<double>("amplitude")*millibarn;
00032 suppressionfactor = p.getParameter<double>("reggeSuppression");
00033
00034 edm::LogInfo("CustomPhysics")<<"Read in physics parameters:"<<G4endl;
00035 edm::LogInfo("CustomPhysics")<<"Resonant = "<< resonant <<G4endl;
00036 edm::LogInfo("CustomPhysics")<<"ResonanceEnergy = "<<ek_0/GeV<<" GeV"<<G4endl;
00037 edm::LogInfo("CustomPhysics")<<"Gamma = "<<gamma/GeV<<" GeV"<<G4endl;
00038 edm::LogInfo("CustomPhysics")<<"Amplitude = "<<amplitude/millibarn<<" millibarn"<<G4endl;
00039 edm::LogInfo("CustomPhysics")<<"ReggeSuppression = "<<100*suppressionfactor<<" %"<<G4endl;
00040
00041 checkfraction = 0;
00042 n_22 = 0;
00043 n_23 = 0;
00044
00045
00046 while(getline(process_stream,line)){
00047 std::vector<G4String> tokens;
00048
00049 ReadAndParse(line,tokens,"#");
00050
00051 G4String incident = tokens[0];
00052
00053 G4ParticleDefinition* incidentDef = particleTable->FindParticle(incident);
00054
00055 G4int incidentPDG = incidentDef->GetPDGEncoding();
00056 known_particles[incidentDef]=true;
00057
00058 G4String target = tokens[1];
00059 edm::LogInfo("CustomPhysics")<<"Incident: "<<incident
00060 <<" Target: "<<target<<G4endl;
00061
00062
00063 ReactionProduct prod;
00064 for (G4int i = 2; i != tokens.size();i++){
00065 G4String part = tokens[i];
00066 if (particleTable->contains(part))
00067 {
00068 prod.push_back(particleTable->FindParticle(part)->GetPDGEncoding());
00069 } else {
00070 edm::LogInfo("CustomPhysics")<<"Particle: "<<part<<" is unknown."<<G4endl;
00071 G4Exception("G4ProcessHelper", "UnkownParticle", FatalException,
00072 "Initialization: The reaction product list contained an unknown particle");
00073 }
00074 }
00075 if (target == "proton")
00076 {
00077 pReactionMap[incidentPDG].push_back(prod);
00078 } else if (target == "neutron") {
00079 nReactionMap[incidentPDG].push_back(prod);
00080 } else {
00081 G4Exception("G4ProcessHelper", "IllegalTarget", FatalException,
00082 "Initialization: The reaction product list contained an illegal target particle");
00083 }
00084 }
00085
00086 process_stream.close();
00087
00088 }
00089
00090 G4bool G4ProcessHelper::ApplicabilityTester(const G4ParticleDefinition& aPart){
00091 const G4ParticleDefinition* aP = &aPart;
00092 if (known_particles[aP]) return true;
00093 return false;
00094 }
00095
00096 G4double G4ProcessHelper::GetInclusiveCrossSection(const G4DynamicParticle *aParticle,
00097 const G4Element *anElement){
00098
00099
00100
00101
00102
00103
00104 G4int thePDGCode = aParticle->GetDefinition()->GetPDGEncoding();
00105
00106 G4double theXsec = 0;
00107
00108
00109 if(CustomPDGParser::s_isRGlueball(thePDGCode)) {
00110 theXsec = 24 * millibarn;
00111 } else {
00112 std::vector<G4int> nq=CustomPDGParser::s_containedQuarks(thePDGCode);
00113
00114 for (std::vector<G4int>::iterator it = nq.begin();
00115 it != nq.end();
00116 it++)
00117 {
00118
00119 if (*it == 1 || *it == 2) theXsec += 12 * millibarn;
00120 if (*it == 3) theXsec += 6 * millibarn;
00121 }
00122 }
00123
00124
00125
00126 if(resonant)
00127 {
00128 double e_0 = ek_0 + aParticle->GetDefinition()->GetPDGMass();
00129
00130 e_0 = sqrt(aParticle->GetDefinition()->GetPDGMass()*aParticle->GetDefinition()->GetPDGMass()
00131 + theProton->GetPDGMass()*theProton->GetPDGMass()
00132 + 2.*e_0*theProton->GetPDGMass());
00133
00134
00135 double sqrts=sqrt(aParticle->GetDefinition()->GetPDGMass()*aParticle->GetDefinition()->GetPDGMass()
00136 + theProton->GetPDGMass()*theProton->GetPDGMass() + 2*aParticle->GetTotalEnergy()*theProton->GetPDGMass());
00137
00138 double res_result = amplitude*(gamma*gamma/4.)/((sqrts-e_0)*(sqrts-e_0)+(gamma*gamma/4.));
00139
00140 theXsec += res_result;
00141
00142 }
00143
00144
00145
00146
00147 return theXsec * pow(anElement->GetN(),0.7)*1.25;
00148
00149 }
00150
00151 ReactionProduct G4ProcessHelper::GetFinalState(const G4Track& aTrack, G4ParticleDefinition*& aTarget){
00152
00153 const G4DynamicParticle* aDynamicParticle = aTrack.GetDynamicParticle();
00154
00155
00156
00157
00158
00159
00160 G4Material* aMaterial = aTrack.GetMaterial();
00161 const G4ElementVector* theElementVector = aMaterial->GetElementVector() ;
00162 const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
00163
00164 G4double NumberOfProtons=0;
00165 G4double NumberOfNucleons=0;
00166
00167 for ( size_t elm=0 ; elm < aMaterial->GetNumberOfElements() ; elm++ )
00168 {
00169
00170 NumberOfProtons += NbOfAtomsPerVolume[elm]*(*theElementVector)[elm]->GetZ();
00171
00172 NumberOfNucleons += NbOfAtomsPerVolume[elm]*(*theElementVector)[elm]->GetN();
00173 }
00174
00175 if(G4UniformRand()<NumberOfProtons/NumberOfNucleons)
00176 {
00177 theReactionMap = &pReactionMap;
00178 theTarget = theProton;
00179 } else {
00180 theReactionMap = &nReactionMap;
00181 theTarget = theNeutron;
00182 }
00183 aTarget = theTarget;
00184 const G4int theIncidentPDG = aDynamicParticle->GetDefinition()->GetPDGEncoding();
00185
00186
00187 ReactionProductList* aReactionProductList = &((*theReactionMap)[theIncidentPDG]);
00188
00189
00190
00191
00192
00193
00194
00195 G4int N22 = 0;
00196 G4int N23 = 0;
00197
00198
00199 ReactionProductList theReactionProductList;
00200 std::vector<bool> theChargeChangeList;
00201
00202 for (ReactionProductList::iterator prod_it = aReactionProductList->begin();
00203 prod_it != aReactionProductList->end();
00204 prod_it++){
00205 G4int secondaries = prod_it->size();
00206
00207 if(ReactionIsPossible(*prod_it,aDynamicParticle)){
00208
00209 theReactionProductList.push_back(*prod_it);
00210 if (secondaries == 2){
00211 N22++;
00212 } else if (secondaries ==3) {
00213 N23++;
00214 } else {
00215 G4cerr << "ReactionProduct has unsupported number of secondaries: "<<secondaries<<G4endl;
00216 }
00217 }
00218
00219
00220 }
00221
00222
00223 if (theReactionProductList.size()==0) G4Exception("G4ProcessHelper", "NoProcessPossible", FatalException,
00224 "GetFinalState: No process could be selected from the given list.");
00225
00226
00227 G4double p22 = 0.15;
00228 G4double p23 = 1-p22;
00229
00230 std::vector<G4double> Probabilities;
00231 std::vector<G4bool> TwotoThreeFlag;
00232
00233 G4double CumulatedProbability = 0;
00234
00235
00236
00237 for (G4int i = 0; std::abs(i) != theReactionProductList.size(); i++){
00238 if (theReactionProductList[i].size() == 2) {
00239 CumulatedProbability += p22/N22;
00240 TwotoThreeFlag.push_back(false);
00241 } else {
00242 CumulatedProbability += p23/N23;
00243 TwotoThreeFlag.push_back(true);
00244 }
00245 Probabilities.push_back(CumulatedProbability);
00246
00247 }
00248
00249
00250
00251 for (std::vector<G4double>::iterator it = Probabilities.begin();
00252 it != Probabilities.end();
00253 it++)
00254 {
00255 *it /= CumulatedProbability;
00256
00257 }
00258
00259
00260
00261
00262 G4bool selected = false;
00263 G4int tries = 0;
00264
00265
00266
00267 G4int i;
00268 while(!selected && tries < 100){
00269 i=0;
00270 G4double dice = G4UniformRand();
00271
00272 while(dice>Probabilities[i] && std::abs(i)<theReactionProductList.size()){
00273
00274 i++;
00275 }
00276
00277
00278
00279 if(!TwotoThreeFlag[i]) {
00280
00281 selected = true;
00282 } else {
00283
00284 if (PhaseSpace(theReactionProductList[i],aDynamicParticle)>G4UniformRand()) selected = true;
00285
00286 }
00287
00288 if(selected&&particleTable->FindParticle(theReactionProductList[i][0])->GetPDGCharge()!=aDynamicParticle->GetDefinition()->GetPDGCharge())
00289 {
00290
00291
00292
00293
00294
00295
00296 if(G4UniformRand()<suppressionfactor) selected = false;
00297 }
00298 tries++;
00299
00300 }
00301 if(tries>=100) G4cerr<<"Could not select process!!!!"<<G4endl;
00302
00303
00304
00305
00306
00307 if (theReactionProductList[i].size()==2) {
00308 n_22++;
00309 } else {
00310 n_23++;
00311 }
00312
00313 checkfraction = (1.0*n_22)/(n_22+n_23);
00314
00315
00316
00317 return theReactionProductList[i];
00318 }
00319
00320 G4double G4ProcessHelper::ReactionProductMass(const ReactionProduct& aReaction,const G4DynamicParticle* aDynamicParticle){
00321
00322 G4double E_incident = aDynamicParticle->GetTotalEnergy();
00323
00324
00325 G4double m_1 = aDynamicParticle->GetDefinition()->GetPDGMass();
00326 G4double m_2 = theTarget->GetPDGMass();
00327
00328 G4double sqrts = sqrt(m_1*m_1 + m_2*(m_2 + 2 * E_incident));
00329
00330
00331 G4double M_after = 0;
00332 for (ReactionProduct::const_iterator r_it = aReaction.begin(); r_it !=aReaction.end(); r_it++){
00333
00334 M_after += particleTable->FindParticle(*r_it)->GetPDGMass();
00335 }
00336
00337 return sqrts - M_after;
00338 }
00339
00340 G4bool G4ProcessHelper::ReactionIsPossible(const ReactionProduct& aReaction,const G4DynamicParticle* aDynamicParticle){
00341 if (ReactionProductMass(aReaction,aDynamicParticle)>0) return true;
00342 return false;
00343 }
00344
00345 G4double G4ProcessHelper::PhaseSpace(const ReactionProduct& aReaction,const G4DynamicParticle* aDynamicParticle){
00346 G4double qValue = ReactionProductMass(aReaction,aDynamicParticle);
00347 G4double phi = sqrt(1+qValue/(2*0.139*GeV))*pow(qValue/(1.1*GeV),3./2.);
00348 return (phi/(1+phi));
00349 }
00350
00351 void G4ProcessHelper::ReadAndParse(const G4String& str,
00352 std::vector<G4String>& tokens,
00353 const G4String& delimiters)
00354 {
00355
00356 G4String::size_type lastPos = str.find_first_not_of(delimiters, 0);
00357
00358 G4String::size_type pos = str.find_first_of(delimiters, lastPos);
00359
00360 while (G4String::npos != pos || G4String::npos != lastPos)
00361 {
00362
00363 G4String temp = str.substr(lastPos, pos - lastPos);
00364 while(temp.c_str()[0] == ' ') temp.erase(0,1);
00365 while(temp[temp.size()-1] == ' ') temp.erase(temp.size()-1,1);
00366
00367 tokens.push_back(temp);
00368
00369 lastPos = str.find_first_not_of(delimiters, pos);
00370
00371 pos = str.find_first_of(delimiters, lastPos);
00372 }
00373 }