G4HadronicProcess.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class source file
00031 //
00032 // G4HadronicProcess
00033 //
00034 // original by H.P.Wellisch
00035 // J.L. Chuma, TRIUMF, 10-Mar-1997
00036 //
00037 // Modifications:
00038 // 05-Jul-2010 V.Ivanchenko cleanup commented lines
00039 // 20-Jul-2011 M.Kelsey -- null-pointer checks in DumpState()
00040 // 24-Sep-2011 M.Kelsey -- Use envvar G4HADRONIC_RANDOM_FILE to save random
00041 //              engine state before each model call
00042 // 18-Oct-2011 M.Kelsey -- Handle final-state cases in conservation checks.
00043 // 14-Mar-2012 G.Folger -- enhance checks for conservation of energy, etc.
00044 // 28-Jul-2012 M.Maire  -- add function GetTargetDefinition() 
00045 // 14-Sep-2012 Inherit from RestDiscrete, use subtype code (now in ctor) to
00046 //              configure base-class
00047 // 28-Sep-2012 Restore inheritance from G4VDiscreteProcess, remove enable-flag
00048 //              changing, remove warning message from original ctor.
00049 
00050 #include "G4HadronicProcess.hh"
00051 
00052 #include "G4Types.hh"
00053 #include "G4SystemOfUnits.hh"
00054 #include "G4HadProjectile.hh"
00055 #include "G4ElementVector.hh"
00056 #include "G4Track.hh"
00057 #include "G4Step.hh"
00058 #include "G4Element.hh"
00059 #include "G4ParticleChange.hh"
00060 #include "G4TransportationManager.hh"
00061 #include "G4Navigator.hh"
00062 #include "G4ProcessVector.hh"
00063 #include "G4ProcessManager.hh"
00064 #include "G4StableIsotopes.hh"
00065 #include "G4HadTmpUtil.hh"
00066 #include "G4NucleiProperties.hh"
00067 
00068 #include "G4HadronicException.hh"
00069 #include "G4HadronicProcessStore.hh"
00070 
00071 #include <typeinfo>
00072 #include <sstream>
00073 #include <iostream>
00074 
00075 #include <stdlib.h>
00076 
00077 // File-scope variable to capture environment variable at startup
00078 
00079 static const char* G4Hadronic_Random_File = getenv("G4HADRONIC_RANDOM_FILE");
00080 
00082 
00083 G4HadronicProcess::G4HadronicProcess(const G4String& processName,
00084                                      G4ProcessType procType)
00085  : G4VDiscreteProcess(processName, procType)
00086 {
00087   SetProcessSubType(fHadronInelastic);  // Default unless subclass changes
00088   
00089   theTotalResult = new G4ParticleChange();
00090   theTotalResult->SetSecondaryWeightByProcess(true);
00091   theInteraction = 0;
00092   theCrossSectionDataStore = new G4CrossSectionDataStore();
00093   G4HadronicProcessStore::Instance()->Register(this);
00094   aScaleFactor = 1;
00095   xBiasOn = false;
00096   G4HadronicProcess_debug_flag = false;
00097 
00098   GetEnergyMomentumCheckEnvvars();
00099 }
00100 
00102 
00103 G4HadronicProcess::G4HadronicProcess(const G4String& processName,
00104                                      G4HadronicProcessType aHadSubType)
00105  : G4VDiscreteProcess(processName, fHadronic)
00106 {
00107   SetProcessSubType(aHadSubType);
00108 
00109   theTotalResult = new G4ParticleChange();
00110   theTotalResult->SetSecondaryWeightByProcess(true);
00111   theInteraction = 0;
00112   theCrossSectionDataStore = new G4CrossSectionDataStore();
00113   G4HadronicProcessStore::Instance()->Register(this);
00114   aScaleFactor = 1;
00115   xBiasOn = false;
00116   G4HadronicProcess_debug_flag = false;
00117 
00118   GetEnergyMomentumCheckEnvvars();
00119 }
00120 
00121 
00122 G4HadronicProcess::~G4HadronicProcess()
00123 {
00124   G4HadronicProcessStore::Instance()->DeRegister(this);
00125   delete theTotalResult;
00126   delete theCrossSectionDataStore;
00127 }
00128 
00129 void G4HadronicProcess::GetEnergyMomentumCheckEnvvars() {
00130   levelsSetByProcess = false;
00131 
00132   epReportLevel = getenv("G4Hadronic_epReportLevel") ?
00133     strtol(getenv("G4Hadronic_epReportLevel"),0,10) : 0;
00134 
00135   epCheckLevels.first = getenv("G4Hadronic_epCheckRelativeLevel") ?
00136     strtod(getenv("G4Hadronic_epCheckRelativeLevel"),0) : DBL_MAX;
00137 
00138   epCheckLevels.second = getenv("G4Hadronic_epCheckAbsoluteLevel") ?
00139     strtod(getenv("G4Hadronic_epCheckAbsoluteLevel"),0) : DBL_MAX;
00140 }
00141 
00142 void G4HadronicProcess::RegisterMe( G4HadronicInteraction *a )
00143 {
00144   if(!a) { return; }
00145   try{GetManagerPointer()->RegisterMe( a );}
00146   catch(G4HadronicException & aE)
00147   {
00148     G4ExceptionDescription ed;
00149     aE.Report(ed);
00150     ed << "Unrecoverable error in " << GetProcessName()
00151        << " to register " << a->GetModelName() << G4endl;
00152     G4Exception("G4HadronicProcess::RegisterMe", "had001", FatalException,
00153                 ed);
00154   }
00155   G4HadronicProcessStore::Instance()->RegisterInteraction(this, a);
00156 }
00157 
00158 void G4HadronicProcess::PreparePhysicsTable(const G4ParticleDefinition& p)
00159 {
00160   if(getenv("G4HadronicProcess_debug")) {
00161     G4HadronicProcess_debug_flag = true;
00162   }
00163   G4HadronicProcessStore::Instance()->RegisterParticle(this, &p);
00164 }
00165 
00166 void G4HadronicProcess::BuildPhysicsTable(const G4ParticleDefinition& p)
00167 {
00168   try
00169   {
00170     theCrossSectionDataStore->BuildPhysicsTable(p);
00171   }
00172   catch(G4HadronicException aR)
00173   {
00174     G4ExceptionDescription ed;
00175     aR.Report(ed);
00176     ed << " hadronic initialisation fails" << G4endl;
00177     G4Exception("G4HadronicProcess::BuildPhysicsTable", "had000", 
00178                 FatalException,ed);
00179   }
00180   G4HadronicProcessStore::Instance()->PrintInfo(&p);
00181 }
00182 
00183 G4double G4HadronicProcess::
00184 GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
00185 {
00186   try
00187   {
00188     theLastCrossSection = aScaleFactor*
00189       theCrossSectionDataStore->GetCrossSection(aTrack.GetDynamicParticle(),
00190                                                 aTrack.GetMaterial());
00191   }
00192   catch(G4HadronicException aR)
00193   {
00194     G4ExceptionDescription ed;
00195     aR.Report(ed);
00196     DumpState(aTrack,"GetMeanFreePath",ed);
00197     ed << " Cross section is not available" << G4endl;
00198     G4Exception("G4HadronicProcess::GetMeanFreePath", "had002", FatalException,
00199                 ed);
00200   }
00201   G4double res = DBL_MAX;
00202   if( theLastCrossSection > 0.0 ) { res = 1.0/theLastCrossSection; }
00203   return res;
00204 }
00205 
00206 G4VParticleChange*
00207 G4HadronicProcess::PostStepDoIt(const G4Track& aTrack, const G4Step&)
00208 {
00209   // if primary is not Alive then do nothing
00210   theTotalResult->Clear();
00211   theTotalResult->Initialize(aTrack);
00212   theTotalResult->ProposeWeight(aTrack.GetWeight());
00213   if(aTrack.GetTrackStatus() != fAlive) { return theTotalResult; }
00214 
00215   // Find cross section at end of step and check if <= 0
00216   //
00217   const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
00218   G4Material* aMaterial = aTrack.GetMaterial();
00219 
00220   G4Element* anElement = 0;
00221   try
00222   {
00223      anElement = theCrossSectionDataStore->SampleZandA(aParticle,
00224                                                        aMaterial,
00225                                                        targetNucleus);
00226   }
00227   catch(G4HadronicException & aR)
00228   {
00229     G4ExceptionDescription ed;
00230     aR.Report(ed);
00231     DumpState(aTrack,"SampleZandA",ed);
00232     ed << " PostStepDoIt failed on element selection" << G4endl;
00233     G4Exception("G4HadronicProcess::PostStepDoIt", "had003", FatalException,
00234                 ed);
00235   }
00236 
00237   // check only for charged particles
00238   if(aParticle->GetDefinition()->GetPDGCharge() != 0.0) {
00239     if (GetElementCrossSection(aParticle, anElement, aMaterial) <= 0.0) {
00240       // No interaction
00241       return theTotalResult;
00242     }    
00243   }
00244 
00245   // Next check for illegal track status
00246   //
00247   if (aTrack.GetTrackStatus() != fAlive && aTrack.GetTrackStatus() != fSuspend) {
00248     if (aTrack.GetTrackStatus() == fStopAndKill ||
00249         aTrack.GetTrackStatus() == fKillTrackAndSecondaries ||
00250         aTrack.GetTrackStatus() == fPostponeToNextEvent) {
00251       G4ExceptionDescription ed;
00252       ed << "G4HadronicProcess: track in unusable state - "
00253          << aTrack.GetTrackStatus() << G4endl;
00254       ed << "G4HadronicProcess: returning unchanged track " << G4endl;
00255       DumpState(aTrack,"PostStepDoIt",ed);
00256       G4Exception("G4HadronicProcess::PostStepDoIt", "had004", JustWarning, ed);
00257     }
00258     // No warning for fStopButAlive which is a legal status here
00259     return theTotalResult;
00260   }
00261 
00262   // Go on to regular case
00263   //
00264   G4double originalEnergy = aParticle->GetKineticEnergy();
00265   G4double kineticEnergy = originalEnergy;
00266 
00267   // Get kinetic energy per nucleon for ions
00268   if(aParticle->GetParticleDefinition()->GetBaryonNumber() > 1.5)
00269           kineticEnergy/=aParticle->GetParticleDefinition()->GetBaryonNumber();
00270 
00271   try
00272   {
00273     theInteraction =
00274       ChooseHadronicInteraction( kineticEnergy, aMaterial, anElement );
00275   }
00276   catch(G4HadronicException & aE)
00277   {
00278     G4ExceptionDescription ed;
00279     aE.Report(ed);
00280     ed << "Target element "<<anElement->GetName()<<"  Z= "
00281        << targetNucleus.GetZ_asInt() << "  A= "
00282        << targetNucleus.GetA_asInt() << G4endl;
00283     DumpState(aTrack,"ChooseHadronicInteraction",ed);
00284     ed << " No HadronicInteraction found out" << G4endl;
00285     G4Exception("G4HadronicProcess::PostStepDoIt", "had005", FatalException,
00286                 ed);
00287   }
00288 
00289   // Initialize the hadronic projectile from the track
00290   thePro.Initialise(aTrack);
00291   G4HadFinalState* result = 0;
00292   G4int reentryCount = 0;
00293 
00294   do
00295   {
00296     try
00297     {
00298       // Save random engine if requested for debugging
00299       if (G4Hadronic_Random_File) {
00300          CLHEP::HepRandom::saveEngineStatus(G4Hadronic_Random_File);
00301       }
00302       // Call the interaction
00303       result = theInteraction->ApplyYourself( thePro, targetNucleus);
00304       ++reentryCount;
00305     }
00306     catch(G4HadronicException aR)
00307     {
00308       G4ExceptionDescription ed;
00309       aR.Report(ed);
00310       ed << "Call for " << theInteraction->GetModelName() << G4endl;
00311       ed << "Target element "<<anElement->GetName()<<"  Z= "
00312          << targetNucleus.GetZ_asInt()
00313          << "  A= " << targetNucleus.GetA_asInt() << G4endl;
00314       DumpState(aTrack,"ApplyYourself",ed);
00315       ed << " ApplyYourself failed" << G4endl;
00316       G4Exception("G4HadronicProcess::PostStepDoIt", "had006", FatalException,
00317                   ed);
00318     }
00319 
00320     // Check the result for catastrophic energy non-conservation
00321     result = CheckResult(thePro,targetNucleus, result);
00322 
00323     if(reentryCount>100) {
00324       G4ExceptionDescription ed;
00325       ed << "Call for " << theInteraction->GetModelName() << G4endl;
00326       ed << "Target element "<<anElement->GetName()<<"  Z= "
00327          << targetNucleus.GetZ_asInt()
00328          << "  A= " << targetNucleus.GetA_asInt() << G4endl;
00329       DumpState(aTrack,"ApplyYourself",ed);
00330       ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
00331       G4Exception("G4HadronicProcess::PostStepDoIt", "had006", FatalException,
00332                   ed);
00333     }
00334   }
00335   while(!result);
00336 
00337   result->SetTrafoToLab(thePro.GetTrafoToLab());
00338 
00339   ClearNumberOfInteractionLengthLeft();
00340 
00341   FillResult(result, aTrack);
00342 
00343   if (epReportLevel != 0) {
00344     CheckEnergyMomentumConservation(aTrack, targetNucleus);
00345   }
00346   return theTotalResult;
00347 }
00348 
00349 
00350 void G4HadronicProcess::ProcessDescription(std::ostream& outFile) const
00351 {
00352   outFile << "The description for this process has not been written yet.\n";
00353 }
00354 
00355 
00356 G4double G4HadronicProcess::XBiasSurvivalProbability()
00357 {
00358   G4double result = 0;
00359   G4double nLTraversed = GetTotalNumberOfInteractionLengthTraversed();
00360   G4double biasedProbability = 1.-std::exp(-nLTraversed);
00361   G4double realProbability = 1-std::exp(-nLTraversed/aScaleFactor);
00362   result = (biasedProbability-realProbability)/biasedProbability;
00363   return result;
00364 }
00365 
00366 G4double G4HadronicProcess::XBiasSecondaryWeight()
00367 {
00368   G4double result = 0;
00369   G4double nLTraversed = GetTotalNumberOfInteractionLengthTraversed();
00370   result =
00371      1./aScaleFactor*std::exp(-nLTraversed/aScaleFactor*(1-1./aScaleFactor));
00372   return result;
00373 }
00374 
00375 void
00376 G4HadronicProcess::FillResult(G4HadFinalState * aR, const G4Track & aT)
00377 {
00378   theTotalResult->ProposeLocalEnergyDeposit(aR->GetLocalEnergyDeposit());
00379 
00380   G4double rotation = CLHEP::twopi*G4UniformRand();
00381   G4ThreeVector it(0., 0., 1.);
00382 
00383   G4double efinal = aR->GetEnergyChange();
00384   if(efinal < 0.0) { efinal = 0.0; }
00385 
00386   // check status of primary
00387   if(aR->GetStatusChange() == stopAndKill) {
00388     theTotalResult->ProposeTrackStatus(fStopAndKill);
00389     theTotalResult->ProposeEnergy( 0.0 );
00390 
00391     // check its final energy
00392   } else if(0.0 == efinal) {
00393     theTotalResult->ProposeEnergy( 0.0 );
00394     if(aT.GetParticleDefinition()->GetProcessManager()
00395        ->GetAtRestProcessVector()->size() > 0)
00396          { aParticleChange.ProposeTrackStatus(fStopButAlive); }
00397     else { aParticleChange.ProposeTrackStatus(fStopAndKill); }
00398 
00399     // primary is not killed apply rotation and Lorentz transformation
00400   } else  {
00401     theTotalResult->ProposeTrackStatus(fAlive);
00402     G4double mass = aT.GetParticleDefinition()->GetPDGMass();
00403     G4double newE = efinal + mass;
00404     G4double newP = std::sqrt(efinal*(efinal + 2*mass));
00405     G4ThreeVector newPV = newP*aR->GetMomentumChange();
00406     G4LorentzVector newP4(newE, newPV);
00407     newP4.rotate(rotation, it);
00408     newP4 *= aR->GetTrafoToLab();
00409     theTotalResult->ProposeMomentumDirection(newP4.vect().unit());
00410     newE = newP4.e() - mass;
00411     if(G4HadronicProcess_debug_flag && newE <= 0.0) {
00412       G4ExceptionDescription ed;
00413       DumpState(aT,"Primary has zero energy after interaction",ed);
00414       G4Exception("G4HadronicProcess::FillResults", "had011", JustWarning, ed);
00415     }
00416     if(newE < 0.0) { newE = 0.0; }
00417     theTotalResult->ProposeEnergy( newE );
00418   }
00419 
00420   // check secondaries: apply rotation and Lorentz transformation
00421   G4int nSec = aR->GetNumberOfSecondaries();
00422   theTotalResult->SetNumberOfSecondaries(nSec);
00423   G4double weight = aT.GetWeight();
00424 
00425   if (nSec > 0) {
00426     G4double time0 = aT.GetGlobalTime();
00427     for (G4int i = 0; i < nSec; ++i) {
00428       G4LorentzVector theM = aR->GetSecondary(i)->GetParticle()->Get4Momentum();
00429       theM.rotate(rotation, it);
00430       theM *= aR->GetTrafoToLab();
00431       aR->GetSecondary(i)->GetParticle()->Set4Momentum(theM);
00432 
00433       // time of interaction starts from zero
00434       G4double time = aR->GetSecondary(i)->GetTime();
00435       if (time < 0.0) { time = 0.0; }
00436 
00437       // take into account global time
00438       time += time0;
00439 
00440       G4Track* track = new G4Track(aR->GetSecondary(i)->GetParticle(),
00441                                    time, aT.GetPosition());
00442       G4double newWeight = weight*aR->GetSecondary(i)->GetWeight();
00443         // G4cout << "#### ParticleDebug "
00444         // <<GetProcessName()<<" "
00445         // <<aR->GetSecondary(i)->GetParticle()->GetDefinition()->GetParticleName()<<" "
00446         // <<aScaleFactor<<" "
00447         // <<XBiasSurvivalProbability()<<" "
00448         // <<XBiasSecondaryWeight()<<" "
00449         // <<aT.GetWeight()<<" "
00450         // <<aR->GetSecondary(i)->GetWeight()<<" "
00451         // <<aR->GetSecondary(i)->GetParticle()->Get4Momentum()<<" "
00452         // <<G4endl;
00453       track->SetWeight(newWeight);
00454       track->SetTouchableHandle(aT.GetTouchableHandle());
00455       theTotalResult->AddSecondary(track);
00456       if (G4HadronicProcess_debug_flag) {
00457         G4double e = track->GetKineticEnergy();
00458         if (e <= 0.0) {
00459           G4ExceptionDescription ed;
00460           DumpState(aT,"Secondary has zero energy",ed);
00461           ed << "Secondary " << track->GetDefinition()->GetParticleName()
00462              << G4endl;
00463           G4Exception("G4HadronicProcess::FillResults", "had011", JustWarning,ed);
00464         }
00465       }
00466     }
00467   }
00468 
00469   aR->Clear();
00470   return;
00471 }
00472 /*
00473 void
00474 G4HadronicProcess::FillTotalResult(G4HadFinalState* aR, const G4Track& aT)
00475 {
00476   theTotalResult->Clear();
00477   theTotalResult->ProposeLocalEnergyDeposit(0.);
00478   theTotalResult->Initialize(aT);
00479   theTotalResult->SetSecondaryWeightByProcess(true);
00480   theTotalResult->ProposeTrackStatus(fAlive);
00481   G4double rotation = CLHEP::twopi*G4UniformRand();
00482   G4ThreeVector it(0., 0., 1.);
00483 
00484   if(aR->GetStatusChange()==stopAndKill)
00485   {
00486     if( xBiasOn && G4UniformRand()<XBiasSurvivalProbability() )
00487     {
00488       theTotalResult->ProposeParentWeight( XBiasSurvivalProbability()*aT.GetWeight() );
00489     }
00490     else
00491     {
00492       theTotalResult->ProposeTrackStatus(fStopAndKill);
00493       theTotalResult->ProposeEnergy( 0.0 );
00494     }
00495   }
00496   else if(aR->GetStatusChange()!=stopAndKill )
00497   {
00498     if(aR->GetStatusChange()==suspend)
00499     {
00500       theTotalResult->ProposeTrackStatus(fSuspend);
00501       if(xBiasOn)
00502       {
00503         G4ExceptionDescription ed;
00504         DumpState(aT,"FillTotalResult",ed);
00505         G4Exception("G4HadronicProcess::FillTotalResult", "had007", FatalException,
00506                     ed,"Cannot cross-section bias a process that suspends tracks.");
00507       }
00508     } else if (aT.GetKineticEnergy() == 0) {
00509       theTotalResult->ProposeTrackStatus(fStopButAlive);
00510     }
00511 
00512     if(xBiasOn && G4UniformRand()<XBiasSurvivalProbability())
00513     {
00514       theTotalResult->ProposeParentWeight( XBiasSurvivalProbability()*aT.GetWeight() );
00515       G4double newWeight = aR->GetWeightChange()*aT.GetWeight();
00516       G4double newM=aT.GetParticleDefinition()->GetPDGMass();
00517       G4double newE=aR->GetEnergyChange() + newM;
00518       G4double newP=std::sqrt(newE*newE - newM*newM);
00519       G4DynamicParticle * aNew =
00520       new G4DynamicParticle(aT.GetParticleDefinition(), newE, newP*aR->GetMomentumChange());
00521       aR->AddSecondary(G4HadSecondary(aNew, newWeight));
00522     }
00523     else
00524     {
00525       G4double newWeight = aR->GetWeightChange()*aT.GetWeight();
00526       theTotalResult->ProposeParentWeight(newWeight); // This is multiplicative
00527       if(aR->GetEnergyChange()>-.5)
00528       {
00529         theTotalResult->ProposeEnergy(aR->GetEnergyChange());
00530       }
00531       G4LorentzVector newDirection(aR->GetMomentumChange().unit(), 1.);
00532       newDirection*=aR->GetTrafoToLab();
00533       theTotalResult->ProposeMomentumDirection(newDirection.vect());
00534     }
00535   }
00536   else
00537   {
00538     G4ExceptionDescription ed;
00539     ed << "Call for " << theInteraction->GetModelName() << G4endl;
00540     ed << "Target Z= " 
00541            << targetNucleus.GetZ_asInt() 
00542            << "  A= " << targetNucleus.GetA_asInt() << G4endl;
00543     DumpState(aT,"FillTotalResult",ed);
00544     G4Exception("G4HadronicProcess", "had008", FatalException,
00545     "use of unsupported track-status.");
00546   }
00547 
00548   if(GetProcessName() != "hElastic" && GetProcessName() != "HadronElastic"
00549      &&  theTotalResult->GetTrackStatus()==fAlive
00550      && aR->GetStatusChange()==isAlive)
00551     {
00552     // Use for debugging:   G4double newWeight = theTotalResult->GetParentWeight();
00553 
00554     G4double newKE = std::max(DBL_MIN, aR->GetEnergyChange());
00555     G4DynamicParticle* aNew = new G4DynamicParticle(aT.GetParticleDefinition(),
00556                                                     aR->GetMomentumChange(),
00557                                                     newKE);
00558     aR->AddSecondary(aNew);
00559     aR->SetStatusChange(stopAndKill);
00560 
00561     theTotalResult->ProposeTrackStatus(fStopAndKill);
00562     theTotalResult->ProposeEnergy( 0.0 );
00563 
00564   }
00565   theTotalResult->ProposeLocalEnergyDeposit(aR->GetLocalEnergyDeposit());
00566   theTotalResult->SetNumberOfSecondaries(aR->GetNumberOfSecondaries());
00567 
00568   if(aR->GetStatusChange() != stopAndKill)
00569   {
00570     G4double newM=aT.GetParticleDefinition()->GetPDGMass();
00571     G4double newE=aR->GetEnergyChange() + newM;
00572     G4double newP=std::sqrt(newE*newE - newM*newM);
00573     G4ThreeVector newPV = newP*aR->GetMomentumChange();
00574     G4LorentzVector newP4(newE, newPV);
00575     newP4.rotate(rotation, it);
00576     newP4*=aR->GetTrafoToLab();
00577     theTotalResult->ProposeMomentumDirection(newP4.vect().unit());
00578   }
00579 
00580   for(G4int i=0; i<aR->GetNumberOfSecondaries(); ++i)
00581   {
00582     G4LorentzVector theM = aR->GetSecondary(i)->GetParticle()->Get4Momentum();
00583     theM.rotate(rotation, it);
00584     theM*=aR->GetTrafoToLab();
00585     aR->GetSecondary(i)->GetParticle()->Set4Momentum(theM);
00586     G4double time = aR->GetSecondary(i)->GetTime();
00587     if(time<0) time = aT.GetGlobalTime();
00588 
00589     G4Track* track = new G4Track(aR->GetSecondary(i)->GetParticle(),
00590                                  time,
00591                                  aT.GetPosition());
00592 
00593     G4double newWeight = aT.GetWeight()*aR->GetSecondary(i)->GetWeight();
00594     if(xBiasOn) { newWeight *= XBiasSecondaryWeight(); }
00595     track->SetWeight(newWeight);
00596     track->SetTouchableHandle(aT.GetTouchableHandle());
00597     theTotalResult->AddSecondary(track);
00598   }
00599 
00600   aR->Clear();
00601   return;
00602 }
00603 */
00604 
00605 void G4HadronicProcess::BiasCrossSectionByFactor(G4double aScale)
00606 {
00607   xBiasOn = true;
00608   aScaleFactor = aScale;
00609   G4String it = GetProcessName();
00610   if( (it != "PhotonInelastic") &&
00611       (it != "ElectroNuclear") &&
00612       (it != "PositronNuclear") )
00613     {
00614       G4ExceptionDescription ed;
00615       G4Exception("G4HadronicProcess::BiasCrossSectionByFactor", "had009", FatalException, ed,
00616                   "Cross-section biasing available only for gamma and electro nuclear reactions.");
00617     }
00618   if(aScale<100)
00619     {
00620       G4ExceptionDescription ed;
00621       G4Exception("G4HadronicProcess::BiasCrossSectionByFactor", "had010", JustWarning,ed,
00622                   "Cross-section bias readjusted to be above safe limit. New value is 100");
00623       aScaleFactor = 100.;
00624     }
00625 }
00626 
00627 G4HadFinalState* G4HadronicProcess::CheckResult(const G4HadProjectile & aPro,const G4Nucleus &aNucleus, G4HadFinalState * result) const
00628 {
00629    // check for catastrophic energy non-conservation, to re-sample the interaction
00630 
00631    G4HadronicInteraction * theModel = GetHadronicInteraction();
00632    G4double nuclearMass(0);
00633    if (theModel){
00634 
00635       // Compute final-state total energy
00636       G4double finalE(0.);
00637       G4int nSec = result->GetNumberOfSecondaries();
00638 
00639       nuclearMass = G4NucleiProperties::GetNuclearMass(aNucleus.GetA_asInt(),
00640                                                        aNucleus.GetZ_asInt());
00641       if (result->GetStatusChange() != stopAndKill) {
00642         // Interaction didn't complete, returned "do nothing" state          => reset nucleus
00643         //  or  the primary survived the interaction (e.g. electro-nuclear ) => keep  nucleus
00644          finalE=result->GetLocalEnergyDeposit() +
00645                 aPro.GetDefinition()->GetPDGMass() + result->GetEnergyChange();
00646          if( nSec == 0 ){
00647             // Since there are no secondaries, there is no recoil nucleus.
00648             // To check energy balance we must neglect the initial nucleus too.
00649             nuclearMass=0.0;
00650          }
00651       }
00652       for (G4int i = 0; i < nSec; i++) {
00653          finalE += result->GetSecondary(i)->GetParticle()->GetTotalEnergy();
00654       }
00655       G4double deltaE= nuclearMass +  aPro.GetTotalEnergy() -  finalE;
00656 
00657       std::pair<G4double, G4double> checkLevels = theModel->GetFatalEnergyCheckLevels();        // (relative, absolute)
00658       if (std::abs(deltaE) > checkLevels.second && std::abs(deltaE) > checkLevels.first*aPro.GetKineticEnergy()){
00659          // do not delete result, this is a pointer to a data member;
00660          result=0;
00661          G4ExceptionDescription desc;
00662          desc << "Warning: Bad energy non-conservation detected, will "
00663               << (epReportLevel<0 ? "abort the event" : "re-sample the interaction") << G4endl
00664               << " Process / Model: " <<  GetProcessName()<< " / " << theModel->GetModelName() << G4endl
00665               << " Primary: " << aPro.GetDefinition()->GetParticleName()
00666               << " (" << aPro.GetDefinition()->GetPDGEncoding() << "),"
00667               << " E= " <<  aPro.Get4Momentum().e()
00668               << ", target nucleus (" << aNucleus.GetZ_asInt() << ","<< aNucleus.GetA_asInt() << ")" << G4endl
00669               << " E(initial - final) = " << deltaE << " MeV." << G4endl;
00670          G4Exception("G4HadronicProcess:CheckResult()", "had012", epReportLevel<0 ? EventMustBeAborted : JustWarning,desc);
00671       }
00672    }
00673    return result;
00674 }
00675 
00676 void
00677 G4HadronicProcess::CheckEnergyMomentumConservation(const G4Track& aTrack,
00678                                                    const G4Nucleus& aNucleus)
00679 {
00680   G4int target_A=aNucleus.GetA_asInt();
00681   G4int target_Z=aNucleus.GetZ_asInt();
00682   G4double targetMass = G4NucleiProperties::GetNuclearMass(target_A,target_Z);
00683   G4LorentzVector target4mom(0, 0, 0, targetMass);
00684 
00685   G4LorentzVector projectile4mom = aTrack.GetDynamicParticle()->Get4Momentum();
00686   G4int track_A = aTrack.GetDefinition()->GetBaryonNumber();
00687   G4int track_Z = G4lrint(aTrack.GetDefinition()->GetPDGCharge());
00688 
00689   G4int initial_A = target_A + track_A;
00690   G4int initial_Z = target_Z + track_Z;
00691 
00692   G4LorentzVector initial4mom = projectile4mom + target4mom;
00693 
00694   // Compute final-state momentum for scattering and "do nothing" results
00695   G4LorentzVector final4mom;
00696   G4int final_A(0), final_Z(0);
00697 
00698   G4int nSec = theTotalResult->GetNumberOfSecondaries();
00699   if (theTotalResult->GetTrackStatus() != fStopAndKill) {  // If it is Alive
00700      // Either interaction didn't complete, returned "do nothing" state
00701      //  or    the primary survived the interaction (e.g. electro-nucleus )
00702      G4Track temp(aTrack);
00703 
00704      // Use the final energy / momentum
00705      temp.SetMomentumDirection(*theTotalResult->GetMomentumDirection());
00706      temp.SetKineticEnergy(theTotalResult->GetEnergy());
00707 
00708      if( nSec == 0 ){
00709         // Interaction didn't complete, returned "do nothing" state
00710         //   - or suppressed recoil  (e.g. Neutron elastic )
00711         final4mom = temp.GetDynamicParticle()->Get4Momentum() + target4mom;
00712         final_A = initial_A;
00713         final_Z = initial_Z;
00714      }else{
00715         // The primary remains in final state (e.g. electro-nucleus )
00716         final4mom = temp.GetDynamicParticle()->Get4Momentum();
00717         final_A = track_A;
00718         final_Z = track_Z;
00719         // Expect that the target nucleus will have interacted,
00720         //  and its products, including recoil, will be included in secondaries.
00721      }
00722   }
00723   if( nSec > 0 ) {
00724     G4Track* sec;
00725 
00726     for (G4int i = 0; i < nSec; i++) {
00727       sec = theTotalResult->GetSecondary(i);
00728       final4mom += sec->GetDynamicParticle()->Get4Momentum();
00729       final_A += sec->GetDefinition()->GetBaryonNumber();
00730       final_Z += G4lrint(sec->GetDefinition()->GetPDGCharge());
00731     }
00732   }
00733 
00734   // Get level-checking information (used to cut-off relative checks)
00735   G4String processName = GetProcessName();
00736   G4HadronicInteraction* theModel = GetHadronicInteraction();
00737   G4String modelName("none");
00738   if (theModel) modelName = theModel->GetModelName();
00739   std::pair<G4double, G4double> checkLevels = epCheckLevels;
00740   if (!levelsSetByProcess) {
00741     if (theModel) checkLevels = theModel->GetEnergyMomentumCheckLevels();
00742     checkLevels.first= std::min(checkLevels.first,  epCheckLevels.first);
00743     checkLevels.second=std::min(checkLevels.second, epCheckLevels.second);
00744   }
00745 
00746   // Compute absolute total-energy difference, and relative kinetic-energy
00747   G4bool checkRelative = (aTrack.GetKineticEnergy() > checkLevels.second);
00748 
00749   G4LorentzVector diff = initial4mom - final4mom;
00750   G4double absolute = diff.e();
00751   G4double relative = checkRelative ? absolute/aTrack.GetKineticEnergy() : 0.;
00752 
00753   G4double absolute_mom = diff.vect().mag();
00754   G4double relative_mom = checkRelative ? absolute_mom/aTrack.GetMomentum().mag() : 0.;
00755 
00756   // Evaluate relative and absolute conservation
00757   G4bool relPass = true;
00758   G4String relResult = "pass";
00759   if (  std::abs(relative) > checkLevels.first
00760          || std::abs(relative_mom) > checkLevels.first) {
00761     relPass = false;
00762     relResult = checkRelative ? "fail" : "N/A";
00763   }
00764 
00765   G4bool absPass = true;
00766   G4String absResult = "pass";
00767   if (   std::abs(absolute) > checkLevels.second
00768       || std::abs(absolute_mom) > checkLevels.second ) {
00769     absPass = false ;
00770     absResult = "fail";
00771   }
00772 
00773   G4bool chargePass = true;
00774   G4String chargeResult = "pass";
00775   if (   (initial_A-final_A)!=0
00776       || (initial_Z-final_Z)!=0 ) {
00777     chargePass = checkLevels.second < DBL_MAX ? false : true;
00778     chargeResult = "fail";
00779    }
00780 
00781   G4bool conservationPass = (relPass || absPass) && chargePass;
00782 
00783   std::stringstream Myout;
00784   G4bool Myout_notempty(false);
00785   // Options for level of reporting detail:
00786   //  0. off
00787   //  1. report only when E/p not conserved
00788   //  2. report regardless of E/p conservation
00789   //  3. report only when E/p not conserved, with model names, process names, and limits
00790   //  4. report regardless of E/p conservation, with model names, process names, and limits
00791   //  negative -1.., as above, but send output to stderr
00792 
00793   if(   std::abs(epReportLevel) == 4
00794         ||      ( std::abs(epReportLevel) == 3 && ! conservationPass ) ){
00795       Myout << " Process: " << processName << " , Model: " <<  modelName << G4endl;
00796       Myout << " Primary: " << aTrack.GetParticleDefinition()->GetParticleName()
00797             << " (" << aTrack.GetParticleDefinition()->GetPDGEncoding() << "),"
00798             << " E= " <<  aTrack.GetDynamicParticle()->Get4Momentum().e()
00799             << ", target nucleus (" << aNucleus.GetZ_asInt() << ","
00800             << aNucleus.GetA_asInt() << ")" << G4endl;
00801       Myout_notempty=true;
00802   }
00803   if (  std::abs(epReportLevel) == 4
00804          || std::abs(epReportLevel) == 2
00805          || ! conservationPass ){
00806 
00807       Myout << "   "<< relResult  <<" relative, limit " << checkLevels.first << ", values E/T(0) = "
00808              << relative << " p/p(0)= " << relative_mom  << G4endl;
00809       Myout << "   "<< absResult << " absolute, limit (MeV) " << checkLevels.second/MeV << ", values E / p (MeV) = "
00810              << absolute/MeV << " / " << absolute_mom/MeV << G4endl;
00811       Myout << "   "<< chargeResult << " charge/baryon number balance " << (initial_Z-final_Z) << " / " << (initial_A-final_A) << " "<<  G4endl;
00812       Myout_notempty=true;
00813 
00814   }
00815   Myout.flush();
00816   if ( Myout_notempty ) {
00817      if (epReportLevel > 0)      G4cout << Myout.str()<< G4endl;
00818      else if (epReportLevel < 0) G4cerr << Myout.str()<< G4endl;
00819   }
00820 }
00821 
00822 
00823 void G4HadronicProcess::DumpState(const G4Track& aTrack,
00824                                   const G4String& method,
00825                                   G4ExceptionDescription& ed)
00826 {
00827   ed << "Unrecoverable error in the method " << method << " of "
00828      << GetProcessName() << G4endl;
00829   ed << "TrackID= "<< aTrack.GetTrackID() << "  ParentID= "
00830      << aTrack.GetParentID()
00831      << "  " << aTrack.GetParticleDefinition()->GetParticleName()
00832      << G4endl;
00833   ed << "Ekin(GeV)= " << aTrack.GetKineticEnergy()/CLHEP::GeV
00834      << ";  direction= " << aTrack.GetMomentumDirection() << G4endl;
00835   ed << "Position(mm)= " << aTrack.GetPosition()/CLHEP::mm << ";";
00836 
00837   if (aTrack.GetMaterial()) {
00838     ed << "  material " << aTrack.GetMaterial()->GetName();
00839   }
00840   ed << G4endl;
00841 
00842   if (aTrack.GetVolume()) {
00843     ed << "PhysicalVolume  <" << aTrack.GetVolume()->GetName()
00844        << ">" << G4endl;
00845   }
00846 }
00847 /* 
00848 G4ParticleDefinition* G4HadronicProcess::GetTargetDefinition()
00849 {
00850   const G4Nucleus* nuc = GetTargetNucleus();
00851   G4int Z = nuc->GetZ_asInt();
00852   G4int A = nuc->GetA_asInt();
00853   return G4ParticleTable::GetParticleTable()->GetIon(Z,A,0*eV);
00854 }
00855 */
00856 /* end of file */

Generated on Mon May 27 17:48:26 2013 for Geant4 by  doxygen 1.4.7