00001
00002 #include "G4ProcessManager.hh"
00003 #include "G4ParticleTable.hh"
00004 #include "G4Track.hh"
00005 #include "G4InelasticInteraction.hh"
00006 #include "Randomize.hh"
00007
00008 #include "SimG4Core/CustomPhysics/interface/ToyModelHadronicProcess.hh"
00009 #include "SimG4Core/CustomPhysics/interface/CustomPDGParser.h"
00010 #include "SimG4Core/CustomPhysics/interface/CustomParticle.h"
00011 #include "SimG4Core/CustomPhysics/interface/HadronicProcessHelper.hh"
00012 #include "SimG4Core/CustomPhysics/interface/Decay3Body.h"
00013
00014 ToyModelHadronicProcess::ToyModelHadronicProcess(HadronicProcessHelper * aHelper, const G4String& processName) :
00015 G4VDiscreteProcess(processName), m_verboseLevel(0), m_helper(aHelper), m_detachCloud(true)
00016 {
00017
00018
00019
00020
00021
00022 }
00023
00024
00025 G4bool ToyModelHadronicProcess::IsApplicable(const G4ParticleDefinition& aP)
00026 {
00027 return m_helper->applicabilityTester(aP);
00028 }
00029
00030 G4double ToyModelHadronicProcess::GetMicroscopicCrossSection(const G4DynamicParticle *particle,
00031 const G4Element *element,
00032 G4double )
00033 {
00034
00035 G4double inclusiveCrossSection = m_helper->inclusiveCrossSection(particle,element);
00036
00037
00038 G4double highestEnergyLimit = 10 * TeV ;
00039 G4double lowestEnergyLimit = 1 * eV;
00040 G4double particleEnergy = particle->GetKineticEnergy();
00041
00042 if (particleEnergy > highestEnergyLimit || particleEnergy < lowestEnergyLimit){
00043 if(m_verboseLevel >= 1) std::cout << "ToyModelHadronicProcess: Energy out of bounds [" <<
00044 lowestEnergyLimit / MeV << "MeV , " <<
00045 highestEnergyLimit / MeV << "MeV ] while it is " << particleEnergy/MeV <<
00046 std::endl;
00047 return 0;
00048 } else {
00049 if(m_verboseLevel >= 3) std::cout << "ToyModelHadronicProcess: Return cross section " << inclusiveCrossSection << std::endl;
00050 return inclusiveCrossSection;
00051 }
00052
00053 }
00054
00055
00056 G4double ToyModelHadronicProcess::GetMeanFreePath(const G4Track& aTrack, G4double, G4ForceCondition*)
00057 {
00058 G4Material *aMaterial = aTrack.GetMaterial();
00059 const G4DynamicParticle *aParticle = aTrack.GetDynamicParticle();
00060 G4double sigma = 0.0;
00061
00062 G4int nElements = aMaterial->GetNumberOfElements();
00063
00064 const G4double *theAtomicNumDensityVector =
00065 aMaterial->GetAtomicNumDensityVector();
00066 G4double aTemp = aMaterial->GetTemperature();
00067
00068 for( G4int i=0; i<nElements; ++i )
00069 {
00070 G4double xSection =
00071 GetMicroscopicCrossSection( aParticle, (*aMaterial->GetElementVector())[i], aTemp);
00072 sigma += theAtomicNumDensityVector[i] * xSection;
00073 }
00074
00075 return 1.0/sigma;
00076
00077 }
00078
00079
00080 G4VParticleChange* ToyModelHadronicProcess::PostStepDoIt(const G4Track& track,
00081 const G4Step& )
00082 {
00083
00084
00085 G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
00086 m_particleChange.Initialize(track);
00087
00088 const G4DynamicParticle* incidentRHadron = track.GetDynamicParticle();
00089 double E_0 = incidentRHadron->GetKineticEnergy();
00090 const G4int theIncidentPDG = incidentRHadron->GetDefinition()->GetPDGEncoding();
00091 const G4ThreeVector aPosition = track.GetPosition();
00092
00093 double gamma = incidentRHadron->GetTotalEnergy()/incidentRHadron->GetDefinition()->GetPDGMass();
00094
00095 CustomParticle* CustomIncident = dynamic_cast<CustomParticle*>(incidentRHadron->GetDefinition());
00096 G4DynamicParticle* cloudParticle = new G4DynamicParticle();
00097
00098
00099
00100
00101 if(CustomIncident==0)
00102 {
00103 std::cout << "ToyModelHadronicProcess::PostStepDoIt Definition of particle cloud not available!!" << std::endl;
00104 } else {
00105 cloudParticle->SetDefinition(CustomIncident->GetCloud());
00106 }
00107
00108
00109 double scale=cloudParticle->GetDefinition()->GetPDGMass()/incidentRHadron->GetDefinition()->GetPDGMass();
00110 G4LorentzVector cloudMomentum;
00111 cloudMomentum.setVectM(incidentRHadron->GetMomentum()*scale,cloudParticle->GetDefinition()->GetPDGMass());
00112 G4LorentzVector gluinoMomentum;
00113 gluinoMomentum.setVectM(incidentRHadron->GetMomentum()*(1.-scale),CustomIncident->GetSpectator()->GetPDGMass());
00114
00115 cloudParticle->Set4Momentum(cloudMomentum);
00116
00117 const G4DynamicParticle *incidentParticle;
00118 if(! m_detachCloud) incidentParticle = incidentRHadron;
00119 else incidentParticle= cloudParticle;
00120
00121 if(m_verboseLevel >= 3)
00122 {
00123 std::cout << "ToyModelHadronicProcess::PostStepDoIt After scaling " << std::endl;
00124 std::cout << " RHadron "<<incidentRHadron->GetDefinition()->GetParticleName()<<" pdgmass= " << incidentRHadron->GetDefinition()->GetPDGMass() / GeV
00125 << " P= " << incidentRHadron->Get4Momentum() / GeV<<" mass= "<< incidentRHadron->Get4Momentum().m() / GeV <<std::endl;
00126 std::cout << " Cloud "<<cloudParticle->GetDefinition()->GetParticleName()<<" pdgmass= " << cloudParticle->GetDefinition()->GetPDGMass() / GeV
00127 << " P= " << cloudParticle->Get4Momentum() / GeV<<" mass= "<< cloudParticle->Get4Momentum().m() / GeV <<std::endl;
00128 std::cout << " Incident pdgmass= " << incidentParticle->GetDefinition()->GetPDGMass() / GeV
00129 << " P= " << incidentParticle->Get4Momentum() / GeV<<" mass= "<< incidentParticle->Get4Momentum().m() / GeV <<std::endl;
00130 }
00131
00132 const G4ThreeVector position = track.GetPosition();
00133 const G4int incidentParticlePDG = incidentRHadron->GetDefinition()->GetPDGEncoding();
00134 std::vector<G4ParticleDefinition*> newParticles;
00135 std::vector<G4ParticleDefinition*> newParticlesRHadron;
00136
00137 G4bool incidentSurvives = false;
00138
00139
00140 G4ParticleDefinition* targetParticle;
00141 HadronicProcessHelper::ReactionProduct reactionProduct = m_helper->finalState(incidentRHadron,track.GetMaterial(),targetParticle);
00142
00143
00144
00145
00146 for(HadronicProcessHelper::ReactionProduct::iterator it = reactionProduct.begin();
00147 it != reactionProduct.end() ;
00148 it++)
00149 {
00150 G4ParticleDefinition * productDefinition =theParticleTable->FindParticle(*it);
00151
00152 if (productDefinition->GetPDGEncoding()==incidentParticlePDG)
00153 incidentSurvives = true;
00154
00155 newParticlesRHadron.push_back(productDefinition);
00156
00157 CustomParticle* cProd = 0;
00158 cProd = dynamic_cast<CustomParticle*>(productDefinition);
00159
00160 if( cProd!=0 && m_detachCloud)
00161 productDefinition=cProd->GetCloud();
00162
00163 newParticles.push_back(productDefinition);
00164 }
00165
00166 int numberOfSecondaries = reactionProduct.size();
00167
00168 if(m_verboseLevel >= 2) std::cout << "ToyModelHadronicProcess::PostStepDoIt N secondaries: "
00169 << numberOfSecondaries << std::endl;
00170
00171
00172
00173
00174
00175
00176
00177 const G4LorentzVector incident4Momentum = incidentParticle->Get4Momentum();
00178 const G4LorentzVector target4Momentum(0,0,0,targetParticle->GetPDGMass());
00179 const G4LorentzVector sum4Momentum = incident4Momentum + target4Momentum;
00180 const G4ThreeVector cmBoost = sum4Momentum.boostVector();
00181 const G4LorentzVector cm4Momentum = sum4Momentum.rest4Vector();
00182
00183 if(m_verboseLevel >= 2) std::cout << "ToyModelHadronicProcess::PostStepDoIt Kinematics in GeV: " << std::endl
00184 << " Boost = " << cmBoost / GeV<< std::endl
00185 << " 4P CM = " << cm4Momentum / GeV << std::endl
00186 << " 4P Inc = " << incident4Momentum / GeV << std::endl
00187 << " 4P Targ = " << target4Momentum / GeV<< std::endl;
00188
00189
00190 const G4double phi_p = 2*pi*G4UniformRand()-pi ;
00191 const G4double theta_p = pi*G4UniformRand() ;
00192 const G4ThreeVector randomDirection(sin(theta_p)*cos(phi_p),
00193 sin(theta_p)*sin(phi_p),
00194 cos(theta_p));
00195
00196
00197 std::vector<G4double> m;
00198 std::vector<G4LorentzVector> fourMomenta;
00199
00200
00201 for(int ip=0;ip<numberOfSecondaries;ip++)
00202 {
00203 m.push_back(newParticles[ip]->GetPDGMass());
00204 }
00205
00206
00207 if (numberOfSecondaries==2){
00208
00209
00210
00211 G4double energy = cm4Momentum.e();
00212 G4ThreeVector p[2];
00213
00214
00215
00216
00217
00218 G4double cmMomentum = 1/(2*energy)*sqrt(energy*energy*energy*energy + m[0]*m[0]*m[0]*m[0] + m[1]*m[1]*m[1]*m[1]
00219 - 2*(m[0]*m[0] + m[1]*m[1])*energy*energy -2*m[0]*m[0]*m[1]*m[1]);
00220 p[0] = cmMomentum * randomDirection;
00221 p[1] = -p[0];
00222
00223 if(m_verboseLevel >= 2) std::cout << "ToyModelHadronicProcess::PostStepDoIt 2->2: " << std::endl
00224 << " Pcm(GeV) = " << cmMomentum / GeV << std::endl;
00225
00226 for(int ip=0;ip<2;ip++)
00227 {
00228
00229 G4double e = sqrt(p[ip].mag2() + m[ip]*m[ip]);
00230
00231 fourMomenta.push_back(G4LorentzVector(p[ip],e));
00232
00233 fourMomenta[ip].boost(cmBoost);
00234
00235 if(m_verboseLevel >= 2) std::cout << " particle " << ip <<" Plab(GeV) = " << fourMomenta[ip] /GeV << std::endl;
00236 }
00237
00238 } else if (numberOfSecondaries==3) {
00239
00240
00241
00242 for (std::vector<G4double>::iterator it=m.begin();it!=m.end();it++)
00243 fourMomenta.push_back(G4LorentzVector(0,0,0,*it));
00244
00245 Decay3Body KinCalc;
00246 KinCalc.doDecay(cm4Momentum, fourMomenta[0], fourMomenta[1], fourMomenta[2] );
00247
00248
00249 HepRotation rotation(randomDirection,G4UniformRand()*2*pi);
00250 for (std::vector<G4LorentzVector>::iterator it = fourMomenta.begin();
00251 it!=fourMomenta.end();
00252 it++)
00253 {
00254 *it *= rotation;
00255 (*it).boost(cmBoost);
00256 }
00257 if(m_verboseLevel >= 3) G4cout<<"Momentum-check: "<<incident4Momentum /GeV<<" GeV vs "
00258 << (fourMomenta[0]+fourMomenta[1]+fourMomenta[2])/GeV<<G4endl;
00259
00260 }
00261
00262
00263
00264 if(incidentSurvives){
00265
00266 m_particleChange.SetNumberOfSecondaries(numberOfSecondaries-1);
00267 if(m_verboseLevel >= 3) std::cout << "Incident survives: set num secondaries to " << numberOfSecondaries-1 << std::endl;
00268
00269 } else {
00270
00271 m_particleChange.SetNumberOfSecondaries(numberOfSecondaries);
00272 m_particleChange.ProposeTrackStatus(fStopAndKill);
00273 if(m_verboseLevel >= 3) std::cout << "Incident does not survive: stopAndKill + set num secondaries to " << numberOfSecondaries << std::endl;
00274 }
00275 double e_kin;
00276
00277 for (int ip=0; ip <numberOfSecondaries;ip++)
00278 {
00279 if (newParticlesRHadron[ip]==incidentRHadron->GetDefinition())
00280 {
00281
00282 if(m_detachCloud)
00283 {
00284 if(m_verboseLevel >= 3)
00285 std::cout << "ToyModelHadronicProcess::PostStepDoIt Add gluino momentum " <<
00286 fourMomenta[ip]/GeV <<"(m="<< fourMomenta[ip].m()/ GeV<<") + " <<gluinoMomentum/GeV<<std::endl; ;
00287 G4LorentzVector p4_new = gluinoMomentum+fourMomenta[ip];
00288 G4ThreeVector momentum = p4_new.vect();
00289 double rhMass=newParticlesRHadron[ip]->GetPDGMass() ;
00290 e_kin = sqrt(momentum.mag()*momentum.mag()+rhMass*rhMass)-rhMass;
00291
00292 fourMomenta[ip].setVectM(momentum,rhMass);
00293
00294 double virt=(p4_new-fourMomenta[ip]).m()/MeV;
00295
00296 if(m_verboseLevel >= 3)
00297 std::cout << " = " << fourMomenta[ip]/GeV <<"(m="<< fourMomenta[ip].m() / GeV<<") vs. "<<rhMass/GeV
00298 << std::endl;
00299 }
00300 m_particleChange.ProposeMomentumDirection(fourMomenta[ip].vect()/fourMomenta[ip].vect().mag());
00301 m_particleChange.ProposeEnergy(fourMomenta[ip].e()-fourMomenta[ip].mag());
00302 if(m_verboseLevel >= 3) std::cout << "ToyModelHadronicProcess::PostStepDoIt Propose momentum " << fourMomenta[ip]/GeV << std::endl;
00303 } else {
00304
00305 G4DynamicParticle* productDynParticle = new G4DynamicParticle();
00306
00307 productDynParticle->SetDefinition(newParticlesRHadron[ip]);
00308
00309 if(newParticlesRHadron[ip]!=newParticles[ip] && m_detachCloud )
00310 {
00311 G4LorentzVector p4_new;
00312 p4_new.setVectM(gluinoMomentum.vect()+fourMomenta[ip].vect(),productDynParticle->GetDefinition()->GetPDGMass());
00313
00314 productDynParticle->Set4Momentum(p4_new);
00315
00316 double virt=(gluinoMomentum+fourMomenta[ip]-p4_new).m()/MeV;
00317
00318 if(m_verboseLevel >= 3)
00319 std::cout << "ToyModelHadronicProcess::PostStepDoIt Add gluino momentum " <<
00320 fourMomenta[ip]/GeV
00321 <<"(m="<< fourMomenta[ip].m()/ GeV<<") + "
00322 <<gluinoMomentum/GeV
00323 <<" = " << productDynParticle->Get4Momentum()/GeV
00324 <<" => "<<productDynParticle->Get4Momentum().m()/GeV
00325 <<" vs. "<<productDynParticle->GetDefinition()->GetPDGMass()/GeV
00326 << std::endl;
00327 }
00328 else
00329 {
00330 productDynParticle->Set4Momentum(fourMomenta[ip]);
00331 }
00332
00333 G4Track* productTrack = new G4Track(productDynParticle,
00334 track.GetGlobalTime(),
00335 position);
00336
00337 if(m_verboseLevel >= 3) std::cout << "ToyModelHadronicProcess::PostStepDoIt Add secondary with 4-Momentum " << productDynParticle->Get4Momentum()/GeV << std::endl;
00338 m_particleChange.AddSecondary(productTrack);
00339 }
00340 }
00341
00342
00343
00344
00345
00346
00347
00348
00349
00350
00351 ClearNumberOfInteractionLengthLeft();
00352
00353 return &m_particleChange;
00354
00355 }
00356
00357 const G4DynamicParticle* ToyModelHadronicProcess::FindRhadron(G4ParticleChange* aParticleChange)
00358 {
00359 G4int nsec = aParticleChange->GetNumberOfSecondaries();
00360 if (nsec==0) return 0;
00361 int i = 0;
00362 G4bool found = false;
00363 while (i!=nsec && !found){
00364 if (dynamic_cast<CustomParticle*>(aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition())!=0) found = true;
00365 i++;
00366 }
00367 i--;
00368 if(found) return aParticleChange->GetSecondary(i)->GetDynamicParticle();
00369 return 0;
00370 }