00001 #include "G4HadReentrentException.hh"
00002 #include "G4ProcessManager.hh"
00003 #include "G4ParticleTable.hh"
00004
00005 #include "SimG4Core/CustomPhysics/interface/FullModelHadronicProcess.hh"
00006 #include "SimG4Core/CustomPhysics/interface/G4ProcessHelper.hh"
00007 #include "SimG4Core/CustomPhysics/interface/Decay3Body.h"
00008 #include "SimG4Core/CustomPhysics/interface/FullModelReactionDynamics.hh"
00009 #include "SimG4Core/CustomPhysics/interface/CustomPDGParser.h"
00010 #include "SimG4Core/CustomPhysics/interface/CustomParticle.h"
00011
00012 FullModelHadronicProcess::FullModelHadronicProcess(G4ProcessHelper * aHelper, const G4String& processName) :
00013 G4VDiscreteProcess(processName), theHelper(aHelper)
00014 {}
00015
00016
00017 FullModelHadronicProcess::~FullModelHadronicProcess(){}
00018
00019 G4bool FullModelHadronicProcess::IsApplicable(const G4ParticleDefinition& aP)
00020 {
00021 return theHelper->ApplicabilityTester(aP);
00022 }
00023
00024 G4double FullModelHadronicProcess::GetMicroscopicCrossSection(const G4DynamicParticle *aParticle,
00025 const G4Element *anElement,
00026 G4double aTemp)
00027 {
00028
00029 G4double InclXsec = theHelper->GetInclusiveCrossSection(aParticle,anElement);
00030
00031
00032
00033 G4double HighestEnergyLimit = 10 * TeV ;
00034 G4double LowestEnergyLimit = 1 * eV;
00035
00036 G4double ParticleEnergy = aParticle->GetKineticEnergy();
00037
00038
00039 if (ParticleEnergy > HighestEnergyLimit || ParticleEnergy < LowestEnergyLimit){
00040 return 0;
00041 } else {
00042
00043 return InclXsec;
00044 }
00045
00046 }
00047
00048 G4double FullModelHadronicProcess::GetMeanFreePath(const G4Track& aTrack, G4double, G4ForceCondition*)
00049 {
00050 G4Material *aMaterial = aTrack.GetMaterial();
00051 const G4DynamicParticle *aParticle = aTrack.GetDynamicParticle();
00052 G4double sigma = 0.0;
00053
00054 G4int nElements = aMaterial->GetNumberOfElements();
00055
00056 const G4double *theAtomicNumDensityVector =
00057 aMaterial->GetAtomicNumDensityVector();
00058 G4double aTemp = aMaterial->GetTemperature();
00059
00060 for( G4int i=0; i<nElements; ++i )
00061 {
00062 G4double xSection =
00063 GetMicroscopicCrossSection( aParticle, (*aMaterial->GetElementVector())[i], aTemp);
00064 sigma += theAtomicNumDensityVector[i] * xSection;
00065 }
00066
00067 return 1.0/sigma;
00068
00069 }
00070
00071 G4VParticleChange* FullModelHadronicProcess::PostStepDoIt(const G4Track& aTrack,
00072 const G4Step& aStep)
00073 {
00074
00075
00076 aParticleChange.Initialize(aTrack);
00077
00078 G4DynamicParticle* IncidentRhadron = const_cast<G4DynamicParticle*>(aTrack.GetDynamicParticle());
00079 CustomParticle* CustomIncident = static_cast<CustomParticle*>(IncidentRhadron->GetDefinition());
00080 const G4ThreeVector aPosition = aTrack.GetPosition();
00081
00082 const G4int theIncidentPDG = IncidentRhadron->GetDefinition()->GetPDGEncoding();
00083 G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
00084 std::vector<G4ParticleDefinition*> theParticleDefinitions;
00085
00086 G4bool IncidentSurvives = false;
00087 G4bool TargetSurvives = false;
00088 G4Nucleus targetNucleus(aTrack.GetMaterial());
00089 G4ParticleDefinition* outgoingRhadron;
00090 G4ParticleDefinition* outgoingCloud;
00091 G4ParticleDefinition* outgoingTarget=0;
00092 double gamma = IncidentRhadron->GetTotalEnergy()/IncidentRhadron->GetDefinition()->GetPDGMass();
00093
00094
00095 G4ThreeVector p_0 = IncidentRhadron->GetMomentum();
00096
00097
00098
00099 G4double e_kin_0 = IncidentRhadron->GetKineticEnergy();
00100
00101
00102 G4DynamicParticle* cloudParticle = new G4DynamicParticle();
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115 cloudParticle->SetDefinition(CustomIncident->GetCloud());
00116
00117 if(cloudParticle->GetDefinition() == 0)
00118 {
00119 std::cout << "ToyModelHadronicProcess::PostStepDoIt Definition of particle cloud not available!!" << std::endl;
00120 }
00121
00122
00123
00124
00125 double scale=cloudParticle->GetDefinition()->GetPDGMass()/IncidentRhadron->GetDefinition()->GetPDGMass();
00126
00127 G4LorentzVector cloudMomentum;
00128 cloudMomentum.setVectM(IncidentRhadron->GetMomentum()*scale,cloudParticle->GetDefinition()->GetPDGMass());
00129 G4LorentzVector gluinoMomentum;
00130
00131 gluinoMomentum.setVectM(IncidentRhadron->GetMomentum()*(1.-scale),CustomIncident->GetSpectator()->GetPDGMass());
00132
00133
00134
00135
00136
00137
00138 G4LorentzVector FullRhadron4Momentum = IncidentRhadron->Get4Momentum();
00139 G4LorentzVector Cloud4Momentum = cloudMomentum;
00140
00141 cloudParticle->Set4Momentum(cloudMomentum);
00142
00143 G4DynamicParticle* OrgPart = cloudParticle;
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153 double E_0 = IncidentRhadron->GetKineticEnergy();
00154
00155
00156 G4double ek = OrgPart->GetKineticEnergy();
00157 G4double amas = OrgPart->GetDefinition()->GetPDGMass();
00158
00159 G4double tkin = targetNucleus.Cinema( ek );
00160 ek += tkin;
00161 OrgPart->SetKineticEnergy( ek );
00162 G4double et = ek + amas;
00163 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00164 G4double pp = OrgPart->GetMomentum().mag();
00165 if( pp > 0.0 )
00166 {
00167 G4ThreeVector momentum = OrgPart->GetMomentum();
00168 OrgPart->SetMomentum( momentum * (p/pp) );
00169 }
00170
00171
00172
00173 tkin = targetNucleus.EvaporationEffects( ek );
00174 ek -= tkin;
00175 if(ek+gluinoMomentum.e()-gluinoMomentum.m()<=0.1*MeV||ek<=0.) {
00176
00177 G4cout<<"Kinetic energy is sick"<<G4endl;
00178 G4cout<<"Full R-hadron: "<<(ek+gluinoMomentum.e()-gluinoMomentum.m())/MeV<<" MeV" <<G4endl;
00179 G4cout<<"Quark system: "<<ek/MeV<<" MeV"<<G4endl;
00180 aParticleChange.ProposeTrackStatus( fStopAndKill );
00181 return &aParticleChange;
00182 }
00183 OrgPart->SetKineticEnergy( ek );
00184 et = ek + amas;
00185 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00186 pp = OrgPart->GetMomentum().mag();
00187
00188 if( pp > 0.0 )
00189 {
00190 G4ThreeVector momentum = OrgPart->GetMomentum();
00191 OrgPart->SetMomentum( momentum * (p/pp) );
00192 }
00193
00194
00195
00196
00197
00198 G4ParticleDefinition* aTarget;
00199 ReactionProduct rp = theHelper->GetFinalState(aTrack,aTarget);
00200 G4bool force2to2 = false;
00201
00202 while(rp.size()!=2 && force2to2){
00203 rp = theHelper->GetFinalState(aTrack,aTarget);
00204 }
00205 G4int NumberOfSecondaries = rp.size();
00206
00207
00208
00209 G4LorentzVector Target4Momentum;
00210 Target4Momentum.setVectM(0.,aTarget->GetPDGMass());
00211
00212 G4LorentzVector psum_full,psum_cloud;
00213 psum_full = FullRhadron4Momentum + Target4Momentum;
00214 psum_cloud = Cloud4Momentum + Target4Momentum;
00215 G4ThreeVector trafo_full_cms = (-1)*psum_full.boostVector();
00216 G4ThreeVector trafo_cloud_cms = (-1)*psum_cloud.boostVector();
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226 for(ReactionProduct::iterator it = rp.begin();
00227 it != rp.end();
00228 it++)
00229 {
00230 G4ParticleDefinition* tempDef = theParticleTable->FindParticle(*it);
00231 CustomParticle* tempCust = dynamic_cast<CustomParticle*>(tempDef);
00232 if (tempDef==aTarget) TargetSurvives = true;
00233
00234
00235 if (tempCust!=0)
00236 {
00237 outgoingRhadron = tempDef;
00238
00239
00240
00241
00242
00243
00244
00245 outgoingCloud=tempCust->GetCloud();
00246 if(outgoingCloud == 0)
00247 {
00248 std::cout << "ToyModelHadronicProcess::PostStepDoIt Definition of outgoing particle cloud not available!!" << std::endl;
00249 }
00250
00251
00252
00253
00254 }
00255
00256 if (tempDef==G4Proton::Proton()||tempDef==G4Neutron::Neutron()) outgoingTarget = tempDef;
00257
00258 if (tempCust==0&&rp.size()==2) outgoingTarget = tempDef;
00259 if (tempDef->GetPDGEncoding()==theIncidentPDG) {
00260 IncidentSurvives = true;
00261 } else {
00262 theParticleDefinitions.push_back(tempDef);
00263
00264
00265
00266
00267
00268 }
00269 }
00270
00271
00272
00273
00274 if (outgoingTarget==0) outgoingTarget = theParticleTable->FindParticle(rp[1]);
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288 if(IncidentSurvives) NumberOfSecondaries--;
00289 aParticleChange.SetNumberOfSecondaries(NumberOfSecondaries);
00290
00291
00292
00293
00294
00295
00296
00297 const G4HadProjectile* originalIncident = new G4HadProjectile(*OrgPart);
00298
00299
00300
00301 G4HadProjectile* org2 = new G4HadProjectile(*OrgPart);
00302 G4LorentzRotation trans = org2->GetTrafoToLab();
00303 delete org2;
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326 G4DynamicParticle *originalTarget = new G4DynamicParticle;
00327 originalTarget->SetDefinition(aTarget);
00328
00329 G4ReactionProduct targetParticle(aTarget);
00330
00331
00332 G4ReactionProduct currentParticle(const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() ) );
00333 currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
00334 currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344 G4ReactionProduct modifiedOriginal = currentParticle;
00345
00346 currentParticle.SetSide( 1 );
00347 targetParticle.SetSide( -1 );
00348 G4bool incidentHasChanged = false;
00349 if (!IncidentSurvives) incidentHasChanged = true;
00350 G4bool targetHasChanged = false;
00351 if (!TargetSurvives) targetHasChanged = true;
00352 G4bool quasiElastic = false;
00353 if (rp.size()==2) quasiElastic = true;
00354 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;
00355 G4int vecLen = 0;
00356 vec.Initialize( 0 );
00357
00358
00359
00360
00361
00362
00363 for (G4int i = 0; i!=NumberOfSecondaries;i++){
00364 if(theParticleDefinitions[i]!=aTarget
00365 && theParticleDefinitions[i]!=originalIncident->GetDefinition()
00366 && theParticleDefinitions[i]!=outgoingRhadron
00367 && theParticleDefinitions[i]!=outgoingTarget)
00368 {
00369 G4ReactionProduct* pa = new G4ReactionProduct;
00370 pa->SetDefinition( theParticleDefinitions[i] );
00371 (G4UniformRand() < 0.5) ? pa->SetSide( -1 ) : pa->SetSide( 1 );
00372 vec.SetElement( vecLen++, pa );
00373 }
00374 }
00375
00376
00377
00378 if (incidentHasChanged) currentParticle.SetDefinitionAndUpdateE( outgoingCloud );
00379 if (incidentHasChanged) modifiedOriginal.SetDefinition( outgoingCloud );
00380 if (targetHasChanged) targetParticle.SetDefinitionAndUpdateE( outgoingTarget );
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390 CalculateMomenta( vec, vecLen,
00391 originalIncident, originalTarget, modifiedOriginal,
00392 targetNucleus, currentParticle, targetParticle,
00393 incidentHasChanged, targetHasChanged, quasiElastic );
00394
00395
00396
00397 G4String cPname = currentParticle.GetDefinition()->GetParticleName();
00398
00399
00400
00401
00402
00403
00404
00405
00406
00407
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426
00427
00428
00429 aParticleChange.SetNumberOfSecondaries(vecLen+NumberOfSecondaries);
00430 G4double e_kin;
00431 G4LorentzVector cloud_p4_new;
00432
00433
00434
00435
00436
00437
00438
00439
00440
00441
00442
00443 cloud_p4_new.setVectM(currentParticle.GetMomentum(),currentParticle.GetMass());
00444 cloud_p4_new *= trans;
00445
00446 G4LorentzVector cloud_p4_old_full = Cloud4Momentum;
00447 cloud_p4_old_full.boost(trafo_full_cms);
00448 G4LorentzVector cloud_p4_old_cloud = Cloud4Momentum;
00449 cloud_p4_old_cloud.boost(trafo_cloud_cms);
00450 G4LorentzVector cloud_p4_full = cloud_p4_new;
00451 cloud_p4_full.boost(trafo_full_cms);
00452 G4LorentzVector cloud_p4_cloud = cloud_p4_new;
00453 cloud_p4_cloud.boost(trafo_cloud_cms);
00454
00455 G4LorentzVector p_g_cms = gluinoMomentum;
00456 p_g_cms.boost(trafo_full_cms);
00457
00458 G4LorentzVector p4_new;
00459
00460
00461
00462 p4_new.setVectM( cloud_p4_new.v() + gluinoMomentum.v(), outgoingRhadron->GetPDGMass() );
00463
00464 double virt=(p4_new-cloud_p4_new-gluinoMomentum).m()/MeV;
00465
00466 G4ThreeVector p_new;
00467 p_new = p4_new.vect();
00468
00469 aParticleChange.ProposeLocalEnergyDeposit((p4_new-cloud_p4_new-gluinoMomentum).m());
00470
00471 if( incidentHasChanged )
00472 {
00473 G4DynamicParticle* p0 = new G4DynamicParticle;
00474 p0->SetDefinition( outgoingRhadron );
00475 p0->SetMomentum( p_new );
00476
00477
00478 G4Track* Track0 = new G4Track(p0,
00479 aTrack.GetGlobalTime(),
00480 aPosition);
00481 aParticleChange.AddSecondary(Track0);
00482
00483
00484
00485
00486
00487 if(p0->GetKineticEnergy()>e_kin_0) {
00488 G4cout<<"ALAAAAARM!!! (incident changed from "
00489 <<IncidentRhadron->GetDefinition()->GetParticleName()
00490 <<" to "<<p0->GetDefinition()->GetParticleName()<<")"<<G4endl;
00491 G4cout<<"Energy loss: "<<(e_kin_0-p0->GetKineticEnergy())/GeV<<" GeV (should be positive)"<<G4endl;
00492
00493 if(rp.size()!=3) G4cout<<"DOUBLE ALAAAAARM!!!"<<G4endl;
00494 }
00495
00496
00497 if(std::abs(p0->GetKineticEnergy()-e_kin_0)>100*GeV) {
00498 G4cout<<"Diff. too big"<<G4endl;
00499 }
00500
00501 aParticleChange.ProposeTrackStatus( fStopAndKill );
00502 }
00503 else
00504 {
00505
00506 G4double p = p_new.mag();
00507 if( p > DBL_MIN )
00508 aParticleChange.ProposeMomentumDirection( p_new.x()/p, p_new.y()/p, p_new.z()/p );
00509 else
00510 aParticleChange.ProposeMomentumDirection( 1.0, 0.0, 0.0 );
00511
00512 G4double aE = sqrt(p*p+(outgoingRhadron->GetPDGMass()*outgoingRhadron->GetPDGMass()) );
00513 e_kin = aE - outgoingRhadron->GetPDGMass();
00514
00515
00516
00517
00518 if(e_kin>e_kin_0) {
00519 G4cout<<"ALAAAAARM!!!"<<G4endl;
00520 G4cout<<"Energy loss: "<<(e_kin_0-e_kin)/GeV<<" GeV (should be positive)"<<G4endl;
00521 if(rp.size()!=3) G4cout<<"DOUBLE ALAAAAARM!!!"<<G4endl;
00522 }
00523 if(std::abs(e_kin-e_kin_0)>100*GeV) {
00524 G4cout<<"Diff. too big"<<G4endl;
00525 }
00526
00527 if (std::fabs(aE)<.1*eV) aE=.1*eV;
00528 aParticleChange.ProposeEnergy( aE- outgoingRhadron->GetPDGMass() );
00529 if(std::abs(e_kin-e_kin_0)>100*GeV) {
00530 G4cout<<"Diff. too big"<<G4endl;
00531 }
00532
00533 }
00534
00535
00536 if( targetParticle.GetMass() > 0.0 )
00537 {
00538 G4DynamicParticle *p1 = new G4DynamicParticle;
00539 p1->SetDefinition( targetParticle.GetDefinition() );
00540
00541 G4ThreeVector momentum = targetParticle.GetMomentum();
00542 momentum = momentum.rotate(cache,what);
00543 p1->SetMomentum( momentum );
00544 p1->SetMomentum( (trans*p1->Get4Momentum()).vect());
00545 G4Track* Track1 = new G4Track(p1,
00546 aTrack.GetGlobalTime(),
00547 aPosition);
00548 aParticleChange.AddSecondary(Track1);
00549 }
00550 G4DynamicParticle *pa;
00551
00552
00553
00554
00555
00556
00557
00558 for(int i=0; i<vecLen; ++i )
00559 {
00560 pa = new G4DynamicParticle();
00561 pa->SetDefinition( vec[i]->GetDefinition() );
00562 pa->SetMomentum( vec[i]->GetMomentum() );
00563 double evec0 = pa->Get4Momentum().e();
00564 pa->Set4Momentum(trans*(pa->Get4Momentum()));
00565 G4ThreeVector pvec = pa->GetMomentum();
00566 double evec = pa->Get4Momentum().e();
00567 G4Track* Trackn = new G4Track(pa,
00568 aTrack.GetGlobalTime(),
00569 aPosition);
00570 aParticleChange.AddSecondary(Trackn);
00571 delete vec[i];
00572 }
00573
00574
00575 const G4DynamicParticle* theRhadron = FindRhadron(&aParticleChange);
00576
00577
00578 if (theRhadron!=NULL||IncidentSurvives)
00579 {
00580
00581 double E_new;
00582 if(IncidentSurvives)
00583 {
00584
00585 E_new = e_kin;
00586 } else {
00587 E_new = theRhadron->GetKineticEnergy();
00588 if(CustomPDGParser::s_isRMeson(theRhadron->GetDefinition()->GetPDGEncoding())
00589 !=CustomPDGParser::s_isRMeson(theIncidentPDG)
00590 ||
00591 CustomPDGParser::s_isMesonino(theRhadron->GetDefinition()->GetPDGEncoding())
00592 !=CustomPDGParser::s_isMesonino(theIncidentPDG)
00593 ) {
00594 std::cout<<"Rm: "<<CustomPDGParser::s_isRMeson(theRhadron->GetDefinition()->GetPDGEncoding())
00595 <<" vs: "<<CustomPDGParser::s_isRMeson(theIncidentPDG)<<std::endl;
00596 std::cout<<"Sm: "<<CustomPDGParser::s_isMesonino(theRhadron->GetDefinition()->GetPDGEncoding())
00597 <<" vs: "<<CustomPDGParser::s_isMesonino(theIncidentPDG)<<std::endl;
00598 }
00599 }
00600
00601
00602 G4LorentzVector p4_old_full = FullRhadron4Momentum;
00603 p4_old_full.boost(trafo_full_cms);
00604 G4LorentzVector p4_old_cloud = FullRhadron4Momentum;
00605 p4_old_cloud.boost(trafo_cloud_cms);
00606 G4LorentzVector p4_full = p4_new;
00607
00608 p4_full=p4_full.boost(trafo_full_cms);
00609
00610 G4LorentzVector p4_cloud = p4_new;
00611 p4_cloud.boost(trafo_cloud_cms);
00612
00613
00614 G4double dtheta_fullcms = p4_full.vect().angle(p4_old_full.vect());
00615 G4double dtheta_cloudcms = p4_cloud.vect().angle(p4_old_cloud.vect());
00616 G4double dtheta_lab = p_new.angle(p_0);
00617
00618 G4double cloud_dtheta_fullcms = cloud_p4_full.vect().angle(cloud_p4_old_full.vect());
00619 G4double cloud_dtheta_cloudcms = cloud_p4_cloud.vect().angle(p4_old_cloud.vect());
00620 G4double cloud_dtheta_lab = cloud_p4_new.vect().angle(p_0);
00621
00622
00623
00624
00625
00626
00627
00628
00629
00630
00631
00632
00633
00634
00635
00636
00637
00638
00639
00640
00641
00642
00643
00644
00645
00646
00647
00648
00649
00650
00651
00652
00653
00654
00655
00656
00657
00658
00659
00660
00661
00662 G4double AbsDeltaE = E_0-E_new;
00663
00664 if (AbsDeltaE > 10*GeV) {
00665 G4cout<<"Energy loss larger than 10 GeV..."<<G4endl;
00666 G4cout<<"E_0: "<<E_0/GeV<<" GeV"<<G4endl;
00667 G4cout<<"E_new: "<<E_new/GeV<<" GeV"<<G4endl;
00668 G4cout<<"Gamma: "<<IncidentRhadron->GetTotalEnergy()/IncidentRhadron->GetDefinition()->GetPDGMass()<<G4endl;
00669 G4cout<<"x: "<<aPosition.x()/cm<<" cm"<<G4endl;
00670 }
00671 }
00672 delete originalIncident;
00673 delete originalTarget;
00674
00675
00676
00677
00678 ClearNumberOfInteractionLengthLeft();
00679
00680
00681 return &aParticleChange;
00682
00683 }
00684
00685
00686 void FullModelHadronicProcess::CalculateMomenta(
00687 G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
00688 G4int &vecLen,
00689 const G4HadProjectile *originalIncident,
00690 const G4DynamicParticle *originalTarget,
00691 G4ReactionProduct &modifiedOriginal,
00692 G4Nucleus &targetNucleus,
00693 G4ReactionProduct ¤tParticle,
00694 G4ReactionProduct &targetParticle,
00695 G4bool &incidentHasChanged,
00696 G4bool &targetHasChanged,
00697 G4bool quasiElastic )
00698 {
00699 FullModelReactionDynamics theReactionDynamics;
00700
00701 cache = 0;
00702 what = originalIncident->Get4Momentum().vect();
00703
00704
00705
00706
00707
00708
00709
00710
00711
00712 if( quasiElastic )
00713 {
00714
00715 theReactionDynamics.TwoBody( vec, vecLen,
00716 modifiedOriginal, originalTarget,
00717 currentParticle, targetParticle,
00718 targetNucleus, targetHasChanged );
00719
00720 return;
00721 }
00722
00723
00724 G4ReactionProduct leadingStrangeParticle;
00725 G4bool leadFlag = MarkLeadingStrangeParticle( currentParticle,
00726 targetParticle,
00727 leadingStrangeParticle );
00728
00729
00730
00731
00732
00733
00734 G4bool finishedGenXPt = false;
00735 G4bool annihilation = false;
00736 if( originalIncident->GetDefinition()->GetPDGEncoding() < 0 &&
00737 currentParticle.GetMass() == 0.0 && targetParticle.GetMass() == 0.0 )
00738 {
00739
00740 annihilation = true;
00741 G4double ekcor = 1.0;
00742 G4double ek = originalIncident->GetKineticEnergy();
00743 G4double ekOrg = ek;
00744
00745 const G4double tarmas = originalTarget->GetDefinition()->GetPDGMass();
00746 if( ek > 1.0*GeV )ekcor = 1./(ek/GeV);
00747 const G4double atomicWeight = targetNucleus.GetN();
00748 ek = 2*tarmas + ek*(1.+ekcor/atomicWeight);
00749
00750 G4double tkin = targetNucleus.Cinema( ek );
00751 ek += tkin;
00752 ekOrg += tkin;
00753 modifiedOriginal.SetKineticEnergy( ekOrg );
00754 }
00755
00756 const G4double twsup[] = { 1.0, 0.7, 0.5, 0.3, 0.2, 0.1 };
00757 G4double rand1 = G4UniformRand();
00758 G4double rand2 = G4UniformRand();
00759 if( annihilation || (vecLen >= 6) ||
00760 (modifiedOriginal.GetKineticEnergy()/GeV >= 1.0) &&
00761 (((originalIncident->GetDefinition() == G4KaonPlus::KaonPlus() ||
00762 originalIncident->GetDefinition() == G4KaonMinus::KaonMinus() ||
00763 originalIncident->GetDefinition() == G4KaonZeroLong::KaonZeroLong() ||
00764 originalIncident->GetDefinition() == G4KaonZeroShort::KaonZeroShort()) &&
00765 rand1 < 0.5) || rand2 > twsup[vecLen]) )
00766 finishedGenXPt =
00767 theReactionDynamics.GenerateXandPt( vec, vecLen,
00768 modifiedOriginal, originalIncident,
00769 currentParticle, targetParticle,
00770 targetNucleus, incidentHasChanged,
00771 targetHasChanged, leadFlag,
00772 leadingStrangeParticle );
00773 if( finishedGenXPt )
00774 {
00775 Rotate(vec, vecLen);
00776 return;
00777 }
00778
00779 G4bool finishedTwoClu = false;
00780 if( modifiedOriginal.GetTotalMomentum()/MeV < 1.0 )
00781 {
00782
00783 for(G4int i=0; i<vecLen; i++) delete vec[i];
00784 vecLen = 0;
00785 }
00786 else
00787 {
00788
00789 theReactionDynamics.SuppressChargedPions( vec, vecLen,
00790 modifiedOriginal, currentParticle,
00791 targetParticle, targetNucleus,
00792 incidentHasChanged, targetHasChanged );
00793
00794 try
00795 {
00796 finishedTwoClu = theReactionDynamics.TwoCluster( vec, vecLen,
00797 modifiedOriginal, originalIncident,
00798 currentParticle, targetParticle,
00799 targetNucleus, incidentHasChanged,
00800 targetHasChanged, leadFlag,
00801 leadingStrangeParticle );
00802 }
00803 catch(G4HadReentrentException aC)
00804 {
00805 aC.Report(G4cout);
00806 throw G4HadReentrentException(__FILE__, __LINE__, "Failing to calculate momenta");
00807 }
00808 }
00809 if( finishedTwoClu )
00810 {
00811 Rotate(vec, vecLen);
00812 return;
00813 }
00814
00815
00816
00817
00818
00819
00820
00821
00822
00823
00824
00825
00826
00827
00828
00829
00830
00831
00832
00833
00834
00835 theReactionDynamics.TwoBody( vec, vecLen,
00836 modifiedOriginal, originalTarget,
00837 currentParticle, targetParticle,
00838 targetNucleus, targetHasChanged );
00839 }
00840
00841
00842 G4bool FullModelHadronicProcess::MarkLeadingStrangeParticle(
00843 const G4ReactionProduct ¤tParticle,
00844 const G4ReactionProduct &targetParticle,
00845 G4ReactionProduct &leadParticle )
00846 {
00847
00848
00849
00850
00851
00852 G4bool lead = false;
00853 if( (currentParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
00854 (currentParticle.GetDefinition() != G4Proton::Proton()) &&
00855 (currentParticle.GetDefinition() != G4Neutron::Neutron()) )
00856 {
00857 lead = true;
00858 leadParticle = currentParticle;
00859 }
00860 else if( (targetParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
00861 (targetParticle.GetDefinition() != G4Proton::Proton()) &&
00862 (targetParticle.GetDefinition() != G4Neutron::Neutron()) )
00863 {
00864 lead = true;
00865 leadParticle = targetParticle;
00866 }
00867 return lead;
00868 }
00869
00870 void FullModelHadronicProcess::Rotate(G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, G4int &vecLen)
00871 {
00872 G4double rotation = 2.*pi*G4UniformRand();
00873 cache = rotation;
00874 G4int i;
00875 for( i=0; i<vecLen; ++i )
00876 {
00877 G4ThreeVector momentum = vec[i]->GetMomentum();
00878 momentum = momentum.rotate(rotation, what);
00879 vec[i]->SetMomentum(momentum);
00880 }
00881 }
00882
00883 const G4DynamicParticle* FullModelHadronicProcess::FindRhadron(G4ParticleChange* aParticleChange)
00884 {
00885 G4int nsec = aParticleChange->GetNumberOfSecondaries();
00886 if (nsec==0) return 0;
00887 int i = 0;
00888 G4bool found = false;
00889 while (i!=nsec && !found){
00890
00891
00892 if (dynamic_cast<CustomParticle*>(aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition())!=0) found = true;
00893 i++;
00894 }
00895 i--;
00896 if(found) return aParticleChange->GetSecondary(i)->GetDynamicParticle();
00897 return 0;
00898 }
00899