CMS 3D CMS Logo

FullModelHadronicProcess.cc

Go to the documentation of this file.
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   //Get the cross section for this particle/element combination from the ProcessHelper
00029   G4double InclXsec = theHelper->GetInclusiveCrossSection(aParticle,anElement);
00030   //  G4cout<<"Returned cross section from helper was: "<<InclXsec/millibarn<<" millibarn"<<G4endl;
00031 
00032   // Need to provide Set-methods for these in time
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     //    G4cout << "Microscopic Cross Section: "<<InclXsec / millibarn<<" millibarn"<<G4endl;
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   //  G4cout<<"*****************    Entering FullModelHadronicProcess::PostStepDoIt       **********************"<<G4endl;
00075   // A little setting up
00076   aParticleChange.Initialize(aTrack);
00077   //  G4DynamicParticle* OrgPart = const_cast<G4DynamicParticle*>(aTrack.GetDynamicParticle());
00078   G4DynamicParticle* IncidentRhadron = const_cast<G4DynamicParticle*>(aTrack.GetDynamicParticle());
00079   CustomParticle* CustomIncident = static_cast<CustomParticle*>(IncidentRhadron->GetDefinition());
00080   const G4ThreeVector aPosition = aTrack.GetPosition();
00081   //  std::cout<<"G: "<<aStep.GetStepLength()/cm<<std::endl;
00082   const G4int theIncidentPDG = IncidentRhadron->GetDefinition()->GetPDGEncoding();
00083   G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
00084   std::vector<G4ParticleDefinition*> theParticleDefinitions;
00085   //  std::vector<G4DynamicParticle*> *theDynamicParticles;//These are probably redundant, but they can easily be removed :-)
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   //  static int n=0;
00098 
00099   G4double e_kin_0 = IncidentRhadron->GetKineticEnergy();
00100   //  G4cout<<e_kin_0/GeV<<G4endl;
00101 
00102   G4DynamicParticle* cloudParticle = new G4DynamicParticle();
00103   /*
00104   if(CustomPDGParser::s_isRMeson(theIncidentPDG))
00105     std::cout<<"Rmeson"<<std::endl;
00106   if(CustomPDGParser::s_isRBaryon(theIncidentPDG)) 
00107     std::cout<<"Rbaryon"<<std::endl;
00108   */
00109   /*
00110   if(CustomPDGParser::s_isRBaryon(theIncidentPDG)) 
00111      cloudParticle->SetDefinition(theParticleTable->FindParticle("rhadronbaryoncloud"));
00112   if(CustomPDGParser::s_isRMeson(theIncidentPDG) || CustomPDGParser::s_isRGlueball(theIncidentPDG) )
00113      cloudParticle->SetDefinition(theParticleTable->FindParticle("rhadronmesoncloud"));
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   std::cout<<"Incoming particle was "<<IncidentRhadron->GetDefinition()->GetParticleName()<<". Corresponding cloud is "<<cloudParticle->GetDefinition()->GetParticleName()<<std::endl;
00123   G4cout<<"Kinetic energy was: "<<IncidentRhadron->GetKineticEnergy()/GeV<<" GeV"<<G4endl;
00124   */
00125   double scale=cloudParticle->GetDefinition()->GetPDGMass()/IncidentRhadron->GetDefinition()->GetPDGMass();
00126   //  std::cout<<"Mass ratio: "<<scale<<std::endl;
00127   G4LorentzVector cloudMomentum;
00128   cloudMomentum.setVectM(IncidentRhadron->GetMomentum()*scale,cloudParticle->GetDefinition()->GetPDGMass());
00129   G4LorentzVector gluinoMomentum;
00130   //  gluinoMomentum.setVectM(IncidentRhadron->GetMomentum()*(1.-scale),theParticleTable->FindParticle("~g")->GetPDGMass()); 
00131   gluinoMomentum.setVectM(IncidentRhadron->GetMomentum()*(1.-scale),CustomIncident->GetSpectator()->GetPDGMass()); 
00132   /*
00133   G4cout <<"Are these the same?"<<G4endl;
00134   G4cout<<gluinoMomentum<<G4endl;
00135   G4cout<<(1.-scale) * IncidentRhadron->Get4Momentum()<<G4endl;
00136   */
00137   //These two for getting CMS transforms later (histogramming purposes...)
00138   G4LorentzVector FullRhadron4Momentum = IncidentRhadron->Get4Momentum();
00139   G4LorentzVector Cloud4Momentum = cloudMomentum;
00140 
00141   cloudParticle->Set4Momentum(cloudMomentum);           
00142 
00143   G4DynamicParticle* OrgPart = cloudParticle;
00144 
00145   /*
00146   std::cout<<"Original momentum: "<<IncidentRhadron->Get4Momentum().v().mag()/GeV<<" GeV, corresponding to gamma: "
00147            <<IncidentRhadron->GetTotalEnergy()/IncidentRhadron->GetDefinition()->GetPDGMass()<<std::endl;
00148   
00149   std::cout<<"Cloud momentum: "<<cloudParticle->Get4Momentum().v().mag()/GeV<<" GeV, corresponding to gamma: "
00150            <<cloudParticle->GetTotalEnergy()/cloudParticle->GetDefinition()->GetPDGMass()<<std::endl;
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   // calculate black track energies
00172   
00173   tkin = targetNucleus.EvaporationEffects( ek );
00174   ek -= tkin;
00175   if(ek+gluinoMomentum.e()-gluinoMomentum.m()<=0.1*MeV||ek<=0.) {
00176     //Very rare event...
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   //Get the final state particles
00197   
00198   G4ParticleDefinition* aTarget; 
00199   ReactionProduct rp = theHelper->GetFinalState(aTrack,aTarget);
00200   G4bool force2to2 = false;
00201   //  G4cout<<"Trying to get final state..."<<G4endl; 
00202   while(rp.size()!=2 && force2to2){
00203     rp = theHelper->GetFinalState(aTrack,aTarget);
00204   }
00205   G4int NumberOfSecondaries = rp.size();
00206   //  G4cout<<"Multiplicity of selected final state: "<<rp.size()<<G4endl;
00207 
00208   //Getting CMS transforms. Boosting is done at histogram filling
00209   G4LorentzVector Target4Momentum;
00210   Target4Momentum.setVectM(0.,aTarget->GetPDGMass());
00211   //  Target4Momentum.setVectM(0.,targetNucleus.GetN()*GeV);
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   psum_full.boost(trafo_full_cms);
00219   psum_cloud.boost(trafo_cloud_cms);
00220   std::cout<<"Checking that the momenta are in deed zero:"<<psum_full.vect()<<std::endl;
00221   */
00222 
00223   // OK Let's make some particles :-)
00224   // We're not using them for anything yet, I know, but let's make sure the machinery is there
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       //      if (tempDef->GetParticleType()=="rhadron")
00235       if (tempCust!=0)
00236         {
00237           outgoingRhadron = tempDef; 
00238           //Setting outgoing cloud definition
00239           /*
00240           if(CustomPDGParser::s_isRBaryon(*it)) 
00241             outgoingCloud=theParticleTable->FindParticle("rhadronbaryoncloud");
00242           if(CustomPDGParser::s_isRMeson(*it) || CustomPDGParser::s_isRGlueball(*it) )
00243             outgoingCloud=theParticleTable->FindParticle("rhadronmesoncloud");
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           std::cout<<"Outgoing Rhadron is: "<<outgoingRhadron->GetParticleName()<<std::endl;
00252           std::cout<<"Outgoing cloud is: "<<outgoingCloud->GetParticleName()<<std::endl;
00253           */
00254         }
00255 
00256       if (tempDef==G4Proton::Proton()||tempDef==G4Neutron::Neutron()) outgoingTarget = tempDef;
00257       //      if (tempDef->GetParticleType()!="rhadron"&&rp.size()==2) outgoingTarget = tempDef;
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         G4DynamicParticle* tempDyn = new G4DynamicParticle();
00265         tempDyn->SetDefinition(tempDef);
00266         theDynamicParticles->push_back(tempDyn);
00267         */
00268       }
00269     }
00270 
00271   //Not using this, so...
00272   //  delete theDynamicParticles;
00273 
00274   if (outgoingTarget==0) outgoingTarget = theParticleTable->FindParticle(rp[1]);
00275 
00276   // A little debug information
00277   /*
00278   G4cout<<"The particles coming out of this reaction will be: ";
00279   for (std::vector<G4DynamicParticle*>::iterator it = theDynamicParticles.begin();
00280        it != theDynamicParticles.end();
00281        it++){
00282     G4cout<< (*it)->GetDefinition()->GetParticleName()<<" ";
00283   }
00284   G4cout<<G4endl;
00285   */
00286   // If the incident particle survives it is not a "secondary", although
00287   // it would probably be easier to fStopAndKill the track every time.
00288   if(IncidentSurvives) NumberOfSecondaries--;
00289   aParticleChange.SetNumberOfSecondaries(NumberOfSecondaries);
00290 
00291 
00292   // ADAPTED FROM G4LEPionMinusInelastic::ApplyYourself
00293   // It is bloody ugly, but such is the way of cut 'n' paste
00294 
00295 
00296   // Set up the incident
00297   const G4HadProjectile* originalIncident = new G4HadProjectile(*OrgPart);//This is where rotation to z-axis is done (Aarrggh!)
00298 
00299 
00300   //Maybe I need to calculate trafo from R-hadron... Bollocks! Labframe does not move. CMS does.  
00301   G4HadProjectile* org2 = new G4HadProjectile(*OrgPart);
00302   G4LorentzRotation trans = org2->GetTrafoToLab();
00303   delete org2;
00304   
00305   //    if (originalIncident->GetKineticEnergy()<= 0.1*MeV) { //Needs rescaling. The kinetic energy of the quarksystem is the relevant quantity
00306   
00307   /*
00308   G4cout<<"Kinetic energies: "<<G4endl;
00309   G4cout<<"True kinetic energy:     "<<originalIncident->GetKineticEnergy()/GeV<<" GeV"<<G4endl;
00310   G4cout<<"Mass:                    "<<originalIncident->GetDefinition()->GetPDGMass()/GeV<<" GeV"<<G4endl;
00311 
00312   G4double e_kin_rescaled = targetNucleus.EvaporationEffects(originalIncident->GetTotalEnergy()-originalIncident->GetDefinition()->GetPDGMass());
00313 
00314   G4cout<<"Rescaled kinetic energy: "<<e_kin_rescaled<<G4endl;
00315 
00316   const G4double cutOff = 0.1*MeV;
00317 
00318   if ( e_kin_rescaled < cutOff )
00319     {
00320       aParticleChange.ProposeTrackStatus( fStopAndKill );//If the dice decides not to cascade I stop the particle
00321       return &aParticleChange;
00322     }
00323   */
00324   // create the target particle
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   G4cout<<"After creation:"<<G4endl;
00338   G4cout<<"currentParticle: "<<currentParticle.GetMomentum()/GeV<<" GeV vs. "<<OrgPart->Get4Momentum()/GeV<<" GeV"<<G4endl;
00339   G4cout<<"targetParticle: "<<targetParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
00340   G4cout<<"Fourmomentum from originalIncident: "<<originalIncident->Get4Momentum()<<" vs "<<OrgPart->Get4Momentum()<<G4endl;
00341   */
00342 
00343 
00344   G4ReactionProduct modifiedOriginal = currentParticle;
00345   
00346   currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
00347   targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
00348   G4bool incidentHasChanged = false;
00349   if (!IncidentSurvives) incidentHasChanged = true; //I wonder if I am supposed to do this...
00350   G4bool targetHasChanged = false;
00351   if (!TargetSurvives) targetHasChanged = true; //Ditto here
00352   G4bool quasiElastic = false;
00353   if (rp.size()==2) quasiElastic = true; //Oh well...
00354   G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;  // vec will contain the secondary particles
00355   G4int vecLen = 0;
00356   vec.Initialize( 0 );
00357 
00358   
00359   
00360   // I hope my understanding of "secondary" is correct here
00361   // I think that it entails only what is not a surviving incident of target
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   //  if (incidentHasChanged) currentParticle.SetDefinitionAndUpdateE( outgoingRhadron );
00377 
00378   if (incidentHasChanged) currentParticle.SetDefinitionAndUpdateE( outgoingCloud );
00379   if (incidentHasChanged) modifiedOriginal.SetDefinition( outgoingCloud );//Is this correct? It solves the "free energy" checking in ReactionDynamics
00380   if (targetHasChanged) targetParticle.SetDefinitionAndUpdateE( outgoingTarget );
00381 
00382   //  G4cout<<"Calling CalculateMomenta... "<<G4endl;
00383   /*
00384   G4cout<<"Current particle starts as: "<<currentParticle.GetDefinition()->GetParticleName()<<G4endl;
00385   G4cout<<"with momentum: "<<currentParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
00386   G4cout<<"May be killed?: "<<currentParticle.GetMayBeKilled()<<G4endl;
00387   */
00388   //  G4double e_temp = currentParticle.GetKineticEnergy();
00389 
00390   CalculateMomenta( vec, vecLen,
00391                     originalIncident, originalTarget, modifiedOriginal,
00392                     targetNucleus, currentParticle, targetParticle,
00393                     incidentHasChanged, targetHasChanged, quasiElastic );
00394 
00395   //  G4cout <<"Cloud loss: "<<(e_temp-currentParticle.GetKineticEnergy())/GeV<<" GeV"<<G4endl;
00396 
00397   G4String cPname = currentParticle.GetDefinition()->GetParticleName();
00398 
00399   //  if(cPname!="rhadronmesoncloud"&&cPname!="rhadronbaryoncloud")
00400   //    {
00401   /*
00402   G4cout<<"Current particle is now: "<<cPname <<G4endl;
00403   G4cout<<"with momentum: "<<currentParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
00404   G4cout<<"and kinetic energy: "<<currentParticle.GetKineticEnergy()/GeV<<" GeV"<<G4endl;
00405   G4cout<<"May be killed?: "<<currentParticle.GetMayBeKilled()<<G4endl;
00406   G4cout<<"Modified original is: "<<modifiedOriginal.GetDefinition()->GetParticleName()<<G4endl;
00407   G4cout<<"with momentum: "<<modifiedOriginal.GetMomentum()/GeV<<" GeV"<<G4endl;
00408   G4cout<<"and kinetic energy: "<<modifiedOriginal.GetKineticEnergy()/GeV<<" GeV"<<G4endl;
00409   G4cout<<"May be killed?: "<<modifiedOriginal.GetMayBeKilled()<<G4endl;
00410   G4cout<<"Target particle is: "<<targetParticle.GetDefinition()->GetParticleName()<<G4endl;
00411   G4cout<<"with momentum: "<<targetParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
00412   G4cout<<"and kinetic energy: "<<targetParticle.GetKineticEnergy()/GeV<<" GeV"<<G4endl;
00413   G4cout<<"May be killed?: "<<targetParticle.GetMayBeKilled()<<G4endl;
00414   G4cout<<"incidentHasChanged: "<<incidentHasChanged<<G4endl;
00415   G4cout<<"targetHasChanged: "<<targetHasChanged<<G4endl;
00416   G4cout<<"Particles in vec:"<<G4endl;
00417   for(int i=0; i<vecLen; ++i )
00418     {
00419       G4cout<< vec[i]->GetDefinition()->GetParticleName()<<G4endl;
00420     }
00421   */
00422 
00423       //    }
00424 
00425   //  G4cout<<"Done!"<<G4endl;
00426   //    const G4LorentzRotation& trans(originalIncident->GetTrafoToLab());
00427   //    G4cout<<"Check aParticleChange.GetNumberOfSecondaries(): "<<aParticleChange.GetNumberOfSecondaries()<<G4endl;
00428 
00429   aParticleChange.SetNumberOfSecondaries(vecLen+NumberOfSecondaries);
00430   G4double e_kin;
00431   G4LorentzVector cloud_p4_new; //Cloud 4-momentum in lab after collision
00432   //  n++;
00433   //  G4cout << n << G4endl;
00434   /*
00435   if(cPname!="rhadronmesoncloud"&&cPname!="rhadronbaryoncloud") {
00436     G4cout<<"Cloud deleted!!! AAARRRRGGGHHH!!!"<<G4endl;
00437     G4cout<<"Cloud name: "<<cPname<<G4endl;
00438     G4cout<<"E_kin_0: "<<e_kin_0/GeV<<" GeV"<<G4endl;
00439     //    G4cout<<"n: "<<n<<G4endl;
00440     //    n=0;
00441   }
00442   */
00443   cloud_p4_new.setVectM(currentParticle.GetMomentum(),currentParticle.GetMass());
00444   cloud_p4_new *= trans;
00445 
00446   G4LorentzVector cloud_p4_old_full = Cloud4Momentum; //quark system in CMS BEFORE collision
00447   cloud_p4_old_full.boost(trafo_full_cms);
00448   G4LorentzVector cloud_p4_old_cloud = Cloud4Momentum; //quark system in cloud CMS BEFORE collision
00449   cloud_p4_old_cloud.boost(trafo_cloud_cms);
00450   G4LorentzVector cloud_p4_full = cloud_p4_new; //quark system in CMS AFTER collision
00451   cloud_p4_full.boost(trafo_full_cms);
00452   G4LorentzVector cloud_p4_cloud = cloud_p4_new; //quark system in cloud CMS AFTER collision
00453   cloud_p4_cloud.boost(trafo_cloud_cms);
00454 
00455   G4LorentzVector p_g_cms = gluinoMomentum; //gluino in CMS BEFORE collision
00456   p_g_cms.boost(trafo_full_cms);
00457 
00458   G4LorentzVector p4_new;
00459   //    p4_new.setVectM(cloud_p4_full.v()+p_g_cms.v(),outgoingRhadron->GetPDGMass());
00460   //    p4_new.boost(-trafo_full_cms);
00461     //    p4_new = cloud_p4_new + gluinoMomentum;
00462   p4_new.setVectM( cloud_p4_new.v() + gluinoMomentum.v(), outgoingRhadron->GetPDGMass() );
00463   //  G4cout<<"P4-diff: "<<(p4_new-cloud_p4_new-gluinoMomentum)/GeV<<", magnitude: "<<(p4_new-cloud_p4_new-gluinoMomentum).m()/MeV<<" MeV" <<G4endl;
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       // May need to run SetDefinitionAndUpdateE here...
00478       G4Track* Track0 = new G4Track(p0,
00479                                     aTrack.GetGlobalTime(),
00480                                     aPosition);
00481       aParticleChange.AddSecondary(Track0);
00482       /*
00483       G4cout<<"Adding a particle "<<p0->GetDefinition()->GetParticleName()<<G4endl;
00484       G4cout<<"with momentum: "<<p0->GetMomentum()/GeV<<" GeV"<<G4endl;
00485       G4cout<<"and kinetic energy: "<<p0->GetKineticEnergy()/GeV<<" GeV"<<G4endl;
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         //Turns out problem is only in 2 -> 3 (Won't fix 2 -> 2 energy deposition)
00493         if(rp.size()!=3) G4cout<<"DOUBLE ALAAAAARM!!!"<<G4endl;
00494       } /*else {
00495         G4cout<<"NO ALAAAAARM!!!"<<G4endl;
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       G4cout<<"New momentum: "<<m/GeV<<" GeV"<<G4endl;
00516       G4cout<<"Kinetic energy: "<<e_kin/GeV<<" GeV"<<G4endl;
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() ); //I do not know if this is correct...
00529       if(std::abs(e_kin-e_kin_0)>100*GeV) {
00530         G4cout<<"Diff. too big"<<G4endl;
00531       }
00532 
00533     }
00534   
00535   //    return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep);
00536   if( targetParticle.GetMass() > 0.0 )  // targetParticle can be eliminated in TwoBody
00537     {
00538       G4DynamicParticle *p1 = new G4DynamicParticle;
00539       p1->SetDefinition( targetParticle.GetDefinition() );
00540       //      G4cout<<"Target secondary: "<<targetParticle.GetDefinition()->GetParticleName()<<G4endl;
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     G4cout<<"vecLen: "<<vecLen<<G4endl;
00553     G4cout<<"#sec's: "<<aParticleChange.GetNumberOfSecondaries()<<G4endl;
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   // Histogram filling  
00575   const G4DynamicParticle* theRhadron = FindRhadron(&aParticleChange);
00576   
00577 
00578   if (theRhadron!=NULL||IncidentSurvives)
00579     {
00580       
00581       double E_new;
00582       if(IncidentSurvives)
00583         {
00584           //      E_new = currentParticle.GetKineticEnergy();
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       //Calculating relevant scattering angles.
00602       G4LorentzVector p4_old_full = FullRhadron4Momentum; //R-hadron in CMS BEFORE collision
00603       p4_old_full.boost(trafo_full_cms);
00604       G4LorentzVector p4_old_cloud = FullRhadron4Momentum; //R-hadron in cloud CMS BEFORE collision
00605       p4_old_cloud.boost(trafo_cloud_cms);
00606       G4LorentzVector p4_full = p4_new; //R-hadron in CMS AFTER collision
00607       //      G4cout<<p4_full.v()/GeV<<G4endl;
00608       p4_full=p4_full.boost(trafo_full_cms);
00609       //      G4cout<<p4_full.m()<<" / "<<(cloud_p4_new+gluinoMomentum).boost(trafo_full_cms).m()<<G4endl;
00610       G4LorentzVector p4_cloud = p4_new; //R-hadron in cloud CMS AFTER collision
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);//acos(p_0*p_new/(p_0.mag()*p_new.mag())); 
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       //Writing out momenta for manual check of boosts:
00623       G4cout<<"******************************************"<<G4endl;
00624 
00625       G4cout<<"R-hadron, before: "<<G4endl;
00626       G4cout<<"Lab: "<<FullRhadron4Momentum.v().mag()/GeV<<" GeV"<<G4endl;
00627       G4cout<<"CMS: "<<p4_old_full.v().mag()/GeV<<" GeV"<<G4endl;
00628       G4cout<<"after: "<<G4endl;
00629       G4cout<<"Lab: "<<p4_new.v().mag()/GeV<<" GeV"<<G4endl;
00630       G4cout<<"CMS: "<<p4_full.v().mag()/GeV<<" GeV"<<G4endl;
00631       G4cout<<"cos(theta*): "<<cos(dtheta_fullcms)<<G4endl;
00632       G4cout<<"Gluino: "<<G4endl;
00633       G4cout<<"Lab: "<<gluinoMomentum.v().mag()/GeV<<" GeV"<<G4endl;
00634       G4cout<<"CMS: "<<p_g_cms.v().mag()/GeV<<" GeV"<<G4endl;
00635 
00636       G4cout<<"Cloud, before: "<<G4endl;
00637       G4cout<<"Lab: "<<Cloud4Momentum.v().mag()/GeV<<" GeV"<<G4endl;
00638       G4cout<<"CMS: "<<cloud_p4_old_full.v().mag()/GeV<<" GeV"<<G4endl;
00639       G4cout<<"after: "<<G4endl;
00640       G4cout<<"CMS: "<<cloud_p4_full.v().mag()/GeV<<" GeV"<<G4endl;
00641       G4cout<<"cloud cos(theta*): "<<cos(cloud_dtheta_fullcms)<<G4endl;
00642       G4cout<<"Longitudinal: "<<cos(cloud_dtheta_fullcms)*cloud_p4_full.v().mag()/GeV<<" GeV"<<G4endl;
00643       G4cout<<"~Combined longitudinal: "<<(cos(cloud_dtheta_fullcms)*cloud_p4_full.v().mag()+p4_old_full.v().mag())/GeV<<" GeV"<<G4endl;
00644       if ((cos(cloud_dtheta_fullcms)*cloud_p4_full.v().mag()+p4_old_full.v().mag())<0.) G4cout<<"ALAAARM"<<G4endl;
00645       G4cout<<"Psum cos(theta*): "<<cos((p_g_cms.v()+cloud_p4_full.v()).angle(p4_old_full.v()))<<G4endl;
00646       G4cout<<"True R-hadron (CMS): "<<(p_g_cms.v()+cloud_p4_full.v()).mag()/GeV<<" GeV"<<G4endl;
00647       G4cout<<"******************************************"<<G4endl;
00648       G4cout<<"Fucking manual fucking calculation:"<<G4endl;
00649       G4cout<<"Momenta:"<<G4endl;
00650       G4cout<<"Cloud, lab: "<<cloud_p4_new.v()/GeV<<" + gluino: "<<gluinoMomentum.v()/GeV
00651             <<" = "<<(cloud_p4_new.v()+gluinoMomentum.v())/GeV
00652             <<". Boost to CMS: "<<(cloud_p4_new+gluinoMomentum).boost(trafo_full_cms).v()/GeV<<G4endl;
00653       G4cout<<"Cloud, CMS: "<<cloud_p4_full.v()/GeV<<" + gluino: "<<p_g_cms.v()/GeV
00654             <<" = "<<(cloud_p4_full.v()+p_g_cms.v())/GeV
00655             <<". Boost to Lab: "<<(cloud_p4_full+p_g_cms).boost(-trafo_full_cms).v()/GeV<<G4endl;
00656       G4cout<<"Ref: "<<p4_new.v()/GeV<<" / "<<p4_full.v()/GeV<<G4endl;
00657       G4cout<<"******************************************"<<G4endl;
00658       */
00659 
00660 
00661       //      std::cout<<"Lab, fullcms, cloudcms: "<<dtheta_lab<<", "<<dtheta_fullcms<<", "<<dtheta_cloudcms<<std::endl;
00662       G4double AbsDeltaE = E_0-E_new;
00663       //      G4cout <<"Energy loss: "<<AbsDeltaE/GeV<<G4endl;
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   //  aParticleChange.DumpInfo();
00675   //  G4cout << "Exiting FullModelHadronicProcess::PostStepDoIt"<<G4endl;
00676 
00677   //clear interaction length      
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,   // the original incident particle
00690                                                const G4DynamicParticle *originalTarget,
00691                                                G4ReactionProduct &modifiedOriginal,   // Fermi motion and evap. effects included
00692                                                G4Nucleus &targetNucleus,
00693                                                G4ReactionProduct &currentParticle,
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   //Commented out like in casqmesmin.F where CALL STPAIR is commented out
00705   /*
00706   theReactionDynamics.ProduceStrangeParticlePairs( vec, vecLen,
00707                                                    modifiedOriginal, originalTarget,
00708                                                    currentParticle, targetParticle,
00709                                                    incidentHasChanged, targetHasChanged );
00710   */
00711 
00712   if( quasiElastic )
00713     {
00714       //      G4cout<<"We are calling TwoBody..."<<G4endl;
00715       theReactionDynamics.TwoBody( vec, vecLen,
00716                                    modifiedOriginal, originalTarget,
00717                                    currentParticle, targetParticle,
00718                                    targetNucleus, targetHasChanged );
00719 
00720       return;
00721     }
00722 
00723   //If ProduceStrangeParticlePairs is commented out, let's cut this one as well
00724   G4ReactionProduct leadingStrangeParticle;
00725   G4bool leadFlag = MarkLeadingStrangeParticle( currentParticle,
00726                                                 targetParticle,
00727                                                 leadingStrangeParticle );
00728 
00729 
00730 
00731   //
00732   // Note: the number of secondaries can be reduced in GenerateXandPt and TwoCluster
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       // original was an anti-particle and annihilation has taken place
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   // PNBlackTrackEnergy is the kinetic energy available for
00817   //   proton/neutron black track particles [was enp(1) in fortran code]
00818   // DTABlackTrackEnergy is the kinetic energy available for
00819   //   deuteron/triton/alpha particles      [was enp(3) in fortran code]
00820   //const G4double pnCutOff = 0.1;
00821   //const G4double dtaCutOff = 0.1;
00822   //if( (targetNucleus.GetN() >= 1.5)
00823   //    && !(incidentHasChanged || targetHasChanged)
00824   //    && (targetNucleus.GetPNBlackTrackEnergy()/MeV <= pnCutOff)
00825   //    && (targetNucleus.GetDTABlackTrackEnergy()/MeV <= dtaCutOff) )
00826   //{
00827   // the atomic weight of the target nucleus is >= 1.5            AND
00828   //   neither the incident nor the target particles have changed  AND
00829   //     there is no kinetic energy available for either proton/neutron
00830   //     or for deuteron/triton/alpha black track particles
00831   // For diffraction scattering on heavy nuclei use elastic routines instead
00832   //G4cerr << "*** Error in G4InelasticInteraction::CalculateMomenta" << G4endl;
00833   //G4cerr << "*** the elastic scattering would be better here ***" <<G4endl;
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 &currentParticle,
00844                                                             const G4ReactionProduct &targetParticle,
00845                                                             G4ReactionProduct &leadParticle )
00846 {
00847   // the following was in GenerateXandPt and TwoCluster
00848   // add a parameter to the GenerateXandPt function telling it about the strange particle
00849   //
00850   // assumes that the original particle was a strange particle
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;              //  set lead to the incident particle
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;              //   set lead to the target particle
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     //    G4cout<<"Checking "<<aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition()->GetParticleName()<<G4endl;
00891     //    if (aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition()->GetParticleType()=="rhadron") found = true;
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 

Generated on Tue Jun 9 17:47:00 2009 for CMSSW by  doxygen 1.5.4