source: trunk/source/processes/hadronic/management/src/G4HadronicProcess.cc@ 1347

Last change on this file since 1347 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 22.0 KB
RevLine 
[819]1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
[1340]26// $Id: G4HadronicProcess.cc,v 1.92 2010/08/17 09:48:20 vnivanch Exp $
27// GEANT4 tag $Name: hadr-man-V09-03-04 $
[819]28//
[1340]29// -------------------------------------------------------------------
30//
31// GEANT4 Class source file
32//
33// G4HadronicProcess
34//
35// original by H.P.Wellisch
36// J.L. Chuma, TRIUMF, 10-Mar-1997
37//
38// Modifications:
39// 05-Jul-2010 V.Ivanchenko cleanup commented lines
40//
[819]41
42#include "G4Types.hh"
43#include "G4HadronicProcess.hh"
[962]44
[819]45#include "G4HadProjectile.hh"
46#include "G4ElementVector.hh"
47#include "G4Track.hh"
48#include "G4Step.hh"
49#include "G4Element.hh"
50#include "G4ParticleChange.hh"
51#include "G4TransportationManager.hh"
52#include "G4Navigator.hh"
53#include "G4ProcessVector.hh"
54#include "G4ProcessManager.hh"
55#include "G4StableIsotopes.hh"
56#include "G4HadTmpUtil.hh"
[1315]57#include "G4NucleiProperties.hh"
[819]58
59#include "G4HadronicException.hh"
[1315]60#include "G4HadronicProcessStore.hh"
[819]61
[1340]62#include <typeinfo>
[962]63
[819]64G4IsoParticleChange * G4HadronicProcess::theIsoResult = 0;
65G4IsoParticleChange * G4HadronicProcess::theOldIsoResult = 0;
66G4bool G4HadronicProcess::isoIsEnabled = true;
67
68void G4HadronicProcess::
69EnableIsotopeProductionGlobally() {isoIsEnabled = true;}
70
71void G4HadronicProcess::
72DisableIsotopeProductionGlobally() {isoIsEnabled = false;}
73
[962]74//////////////////////////////////////////////////////////////////
75
[1315]76G4HadronicProcess::G4HadronicProcess(const G4String& processName,
77 G4ProcessType aType)
78 :G4VDiscreteProcess(processName, aType)
79{
[819]80 ModelingState = 0;
81 isoIsOnAnyway = -1;
82 theTotalResult = new G4ParticleChange();
[1315]83 theInteraction = 0;
[819]84 theCrossSectionDataStore = new G4CrossSectionDataStore();
[962]85 G4HadronicProcessStore::Instance()->Register(this);
[819]86 aScaleFactor = 1;
87 xBiasOn = false;
[962]88 G4HadronicProcess_debug_flag = false;
[1315]89 epReportLevel = 0;
90 epCheckLevels.first = DBL_MAX;
91 epCheckLevels.second = DBL_MAX;
92 levelsSetByProcess = false;
[819]93}
94
95G4HadronicProcess::~G4HadronicProcess()
96{
[962]97 G4HadronicProcessStore::Instance()->DeRegister(this);
[819]98 delete theTotalResult;
99
100 std::for_each(theProductionModels.begin(),
101 theProductionModels.end(), G4Delete());
102
[962]103 delete theOldIsoResult;
104 delete theIsoResult;
[819]105 delete theCrossSectionDataStore;
106}
107
108void G4HadronicProcess::RegisterMe( G4HadronicInteraction *a )
109{
[1315]110 if(!a) { return; }
[819]111 try{GetManagerPointer()->RegisterMe( a );}
112 catch(G4HadronicException & aE)
113 {
[1315]114 G4cout << "Unrecoverable error in " << GetProcessName()
115 << " to register " << a->GetModelName() << G4endl;
[819]116 G4Exception("G4HadronicProcess", "007", FatalException,
117 "Could not register G4HadronicInteraction");
118 }
[962]119 G4HadronicProcessStore::Instance()->RegisterInteraction(this, a);
[819]120}
121
[962]122void G4HadronicProcess::PreparePhysicsTable(const G4ParticleDefinition& p)
123{
[1340]124 if(getenv("G4HadronicProcess_debug")) {
125 G4HadronicProcess_debug_flag = true;
126 }
[962]127 G4HadronicProcessStore::Instance()->RegisterParticle(this, &p);
128}
129
130void G4HadronicProcess::BuildPhysicsTable(const G4ParticleDefinition& p)
131{
132 theCrossSectionDataStore->BuildPhysicsTable(p);
133 G4HadronicProcessStore::Instance()->PrintInfo(&p);
134}
135
[819]136G4double G4HadronicProcess::
137GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
138{
139 try
140 {
[962]141 theLastCrossSection = aScaleFactor*
142 theCrossSectionDataStore->GetCrossSection(aTrack.GetDynamicParticle(),
143 aTrack.GetMaterial());
[819]144 }
145 catch(G4HadronicException aR)
146 {
[1315]147 DumpState(aTrack,"GetMeanFreePath");
[819]148 G4Exception("G4HadronicProcess", "007", FatalException,
149 "G4HadronicProcess::GetMeanFreePath failed");
150 }
[962]151 G4double res = DBL_MAX;
[1315]152 if( theLastCrossSection > 0.0 ) { res = 1.0/theLastCrossSection; }
[962]153 return res;
[819]154}
155
[962]156G4double G4HadronicProcess::
157GetMicroscopicCrossSection(const G4DynamicParticle *aParticle,
158 const G4Element *anElement,
159 G4double aTemp )
[819]160{
[1340]161 G4double x =
[962]162 theCrossSectionDataStore->GetCrossSection(aParticle, anElement, aTemp);
[1340]163 if(x < 0.0) { x = 0.0; }
164 return x;
[819]165}
166
[962]167G4VParticleChange *G4HadronicProcess::PostStepDoIt(
168 const G4Track &aTrack, const G4Step &)
[819]169{
[962]170 if(aTrack.GetTrackStatus() != fAlive && aTrack.GetTrackStatus() != fSuspend) {
171 if (aTrack.GetTrackStatus() == fStopAndKill ||
172 aTrack.GetTrackStatus() == fKillTrackAndSecondaries ||
173 aTrack.GetTrackStatus() == fPostponeToNextEvent) {
[1315]174 G4cout << "G4HadronicProcess: track in unusable state - "
[962]175 << aTrack.GetTrackStatus() << G4endl;
[1315]176 G4cout << "G4HadronicProcess: returning unchanged track " << G4endl;
177 DumpState(aTrack,"PostStepDoIt");
[962]178 G4Exception("G4HadronicProcess", "001", JustWarning, "bailing out");
179 }
180 // No warning for fStopButAlive which is a legal status here
[819]181 theTotalResult->Clear();
182 theTotalResult->Initialize(aTrack);
183 return theTotalResult;
184 }
185
[962]186 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
187 G4Material* aMaterial = aTrack.GetMaterial();
[819]188 G4double originalEnergy = aParticle->GetKineticEnergy();
189 G4double kineticEnergy = originalEnergy;
190
191 // Get kinetic energy per nucleon for ions
[1340]192 if(aParticle->GetParticleDefinition()->GetBaryonNumber() > 1.5)
193 kineticEnergy/=aParticle->GetParticleDefinition()->GetBaryonNumber();
[819]194
195 G4Element* anElement = 0;
196 try
197 {
[962]198 // anElement = ChooseAandZ( aParticle, aMaterial );
199 anElement = theCrossSectionDataStore->SampleZandA(aParticle,
200 aMaterial,
201 targetNucleus);
[819]202 }
203 catch(G4HadronicException & aR)
204 {
[1315]205 DumpState(aTrack,"SampleZandA");
[819]206 G4Exception("G4HadronicProcess", "007", FatalException,
[962]207 "PostStepDoIt failed on element selection.");
[819]208 }
209
210 try
211 {
[1315]212 theInteraction =
213 ChooseHadronicInteraction( kineticEnergy, aMaterial, anElement );
[819]214 }
215 catch(G4HadronicException & aE)
216 {
[1315]217 G4cout << "Target element "<<anElement->GetName()<<" Z= "
218 << targetNucleus.GetZ() << " A= " << targetNucleus.GetN() << G4endl;
219 DumpState(aTrack,"ChooseHadronicInteraction");
[819]220 G4Exception("G4HadronicProcess", "007", FatalException,
221 "ChooseHadronicInteraction failed.");
222 }
223
224 // Initialize the hadronic projectile from the track
225
226 G4HadProjectile thePro(aTrack);
227
228 G4HadFinalState* result = 0;
229 G4int reentryCount = 0;
230
231 do
232 {
233 try
234 {
235 // Call the interaction
[962]236 result = theInteraction->ApplyYourself( thePro, targetNucleus);
[1315]237 ++reentryCount;
[819]238 }
239 catch(G4HadronicException aR)
240 {
[1315]241 G4cout << "Call for " << theInteraction->GetModelName() << G4endl;
242 G4cout << "Target element "<<anElement->GetName()<<" Z= "
243 << targetNucleus.GetZ() << " A= " << targetNucleus.GetN() << G4endl;
244 DumpState(aTrack,"ApplyYourself");
[819]245 G4Exception("G4HadronicProcess", "007", FatalException,
[962]246 "PostStepDoIt failed.");
[819]247 }
[1315]248 if(reentryCount>100) {
249 G4cout << "Call for " << theInteraction->GetModelName() << G4endl;
250 G4cout << "Target element "<<anElement->GetName()<<" Z= "
251 << targetNucleus.GetZ() << " A= " << targetNucleus.GetN() << G4endl;
252 DumpState(aTrack,"ApplyYourself");
253 G4Exception("G4HadronicProcess", "007", FatalException,
254 "Reentering ApplyYourself too often - PostStepDoIt failed.");
255 }
[819]256 }
257 while(!result);
258
259 result->SetTrafoToLab(thePro.GetTrafoToLab());
260
261 ClearNumberOfInteractionLengthLeft();
262 if(isoIsOnAnyway!=-1)
263 {
264 if(isoIsEnabled||isoIsOnAnyway)
265 {
266 result = DoIsotopeCounting(result, aTrack, targetNucleus);
267 }
268 }
269
270 // Put hadronic final state particles into G4ParticleChange
271
272 FillTotalResult(result, aTrack);
[1315]273
[1340]274 if (epReportLevel > 0) {
275 CheckEnergyMomentumConservation(aTrack, targetNucleus);
276 }
[819]277 return theTotalResult;
278}
279
280G4HadFinalState*
281G4HadronicProcess::DoIsotopeCounting(G4HadFinalState * aResult,
282 const G4Track & aTrack,
283 const G4Nucleus & aNucleus)
284{
285 // get the PC from iso-production
286 delete theOldIsoResult;
287 theOldIsoResult = 0;
288 delete theIsoResult;
289 theIsoResult = new G4IsoParticleChange;
290 G4bool done = false;
291 G4IsoResult * anIsoResult = 0;
292 for(unsigned int i=0; i<theProductionModels.size(); i++)
293 {
294 anIsoResult = theProductionModels[i]->GetIsotope(aTrack, aNucleus);
295 if(anIsoResult!=0)
296 {
297 done = true;
298 break;
299 }
300 }
301
302 // If no production models active, use default iso production
303 if(!done) anIsoResult = ExtractResidualNucleus(aTrack, aNucleus, aResult);
304
305 // Add all info explicitely and add typename from model called.
306 theIsoResult->SetIsotope(anIsoResult->GetIsotope());
307 theIsoResult->SetProductionPosition(aTrack.GetPosition());
308 theIsoResult->SetProductionTime(aTrack.GetGlobalTime());
309 theIsoResult->SetParentParticle(*aTrack.GetDynamicParticle());
310 theIsoResult->SetMotherNucleus(anIsoResult->GetMotherNucleus());
311 theIsoResult->SetProducer(typeid(*theInteraction).name());
312
313 delete anIsoResult;
314
315 // If isotope production is enabled the GetIsotopeProductionInfo()
316 // method must be called or else a memory leak will result
317 //
318 // The following code will fix the memory leak, but remove the
319 // isotope information:
320 //
321 // if(theIsoResult) {
322 // delete theIsoResult;
323 // theIsoResult = 0;
324 // }
325
326 return aResult;
327}
328
329G4IsoResult*
330G4HadronicProcess::ExtractResidualNucleus(const G4Track&,
331 const G4Nucleus& aNucleus,
332 G4HadFinalState* aResult)
333{
334 G4double A = aNucleus.GetN();
335 G4double Z = aNucleus.GetZ();
336 G4double bufferA = 0;
337 G4double bufferZ = 0;
338
339 // loop over aResult, and decrement A, Z accordingly
340 // cash the max
[1340]341 for(G4int i=0; i<aResult->GetNumberOfSecondaries(); ++i)
[819]342 {
343 G4HadSecondary* aSecTrack = aResult->GetSecondary(i);
[1340]344 const G4ParticleDefinition* part = aSecTrack->GetParticle()->GetParticleDefinition();
345 G4double Q = part->GetPDGCharge()/eplus;
346 G4double N = part->GetBaryonNumber();
347 if(bufferA < N)
[819]348 {
[1340]349 bufferA = N;
350 bufferZ = Q;
[819]351 }
[1340]352 Z -= Q;
353 A -= N;
[819]354 }
355
356 // if the fragment was part of the final state, it is
357 // assumed to be the heaviest secondary.
358 if(A<0.1)
359 {
360 A = bufferA;
361 Z = bufferZ;
362 }
363
364 // prepare the IsoResult.
365
366 std::ostringstream ost1;
367 ost1 <<Z<<"_"<<A;
368 G4String biff = ost1.str();
369 G4IsoResult * theResult = new G4IsoResult(biff, aNucleus);
370
371 return theResult;
372}
373
374G4double G4HadronicProcess::XBiasSurvivalProbability()
375{
376 G4double result = 0;
377 G4double nLTraversed = GetTotalNumberOfInteractionLengthTraversed();
378 G4double biasedProbability = 1.-std::exp(-nLTraversed);
379 G4double realProbability = 1-std::exp(-nLTraversed/aScaleFactor);
380 result = (biasedProbability-realProbability)/biasedProbability;
381 return result;
382}
383
384G4double G4HadronicProcess::XBiasSecondaryWeight()
385{
386 G4double result = 0;
387 G4double nLTraversed = GetTotalNumberOfInteractionLengthTraversed();
388 result =
389 1./aScaleFactor*std::exp(-nLTraversed/aScaleFactor*(1-1./aScaleFactor));
390 return result;
391}
392
393void
394G4HadronicProcess::FillTotalResult(G4HadFinalState * aR, const G4Track & aT)
395{
[962]396 // G4Nancheck go_wild;
[819]397 theTotalResult->Clear();
398 theTotalResult->ProposeLocalEnergyDeposit(0.);
399 theTotalResult->Initialize(aT);
400 theTotalResult->SetSecondaryWeightByProcess(true);
401 theTotalResult->ProposeTrackStatus(fAlive);
[1315]402 G4double rotation = CLHEP::twopi*G4UniformRand();
[819]403 G4ThreeVector it(0., 0., 1.);
404
405 if(aR->GetStatusChange()==stopAndKill)
406 {
407 if( xBiasOn && G4UniformRand()<XBiasSurvivalProbability() )
408 {
409 theTotalResult->ProposeParentWeight( XBiasSurvivalProbability()*aT.GetWeight() );
410 }
411 else
412 {
413 theTotalResult->ProposeTrackStatus(fStopAndKill);
414 theTotalResult->ProposeEnergy( 0.0 );
415 }
416 }
417 else if(aR->GetStatusChange()!=stopAndKill )
418 {
419 if(aR->GetStatusChange()==suspend)
420 {
421 theTotalResult->ProposeTrackStatus(fSuspend);
422 if(xBiasOn)
423 {
[1315]424 DumpState(aT,"FillTotalResult");
[819]425 G4Exception("G4HadronicProcess", "007", FatalException,
426 "Cannot cross-section bias a process that suspends tracks.");
427 }
428 } else if (aT.GetKineticEnergy() == 0) {
429 theTotalResult->ProposeTrackStatus(fStopButAlive);
430 }
431
432 if(xBiasOn && G4UniformRand()<XBiasSurvivalProbability())
433 {
434 theTotalResult->ProposeParentWeight( XBiasSurvivalProbability()*aT.GetWeight() );
435 G4double newWeight = aR->GetWeightChange()*aT.GetWeight();
[1340]436 G4double newM=aT.GetParticleDefinition()->GetPDGMass();
[819]437 G4double newE=aR->GetEnergyChange() + newM;
438 G4double newP=std::sqrt(newE*newE - newM*newM);
439 G4DynamicParticle * aNew =
[1340]440 new G4DynamicParticle(aT.GetParticleDefinition(), newE, newP*aR->GetMomentumChange());
[819]441 G4HadSecondary * theSec = new G4HadSecondary(aNew, newWeight);
442 aR->AddSecondary(theSec);
443 }
444 else
445 {
446 G4double newWeight = aR->GetWeightChange()*aT.GetWeight();
447 theTotalResult->ProposeParentWeight(newWeight); // This is multiplicative
448 if(aR->GetEnergyChange()>-.5)
449 {
450 theTotalResult->ProposeEnergy(aR->GetEnergyChange());
451 }
452 G4LorentzVector newDirection(aR->GetMomentumChange().unit(), 1.);
453 newDirection*=aR->GetTrafoToLab();
454 theTotalResult->ProposeMomentumDirection(newDirection.vect());
455 }
456 }
457 else
458 {
[1315]459 G4cout << "Call for " << theInteraction->GetModelName() << G4endl;
460 G4cout << "Target Z= "
461 << targetNucleus.GetZ() << " A= " << targetNucleus.GetN() << G4endl;
462 DumpState(aT,"FillTotalResult");
[819]463 G4Exception("G4HadronicProcess", "007", FatalException,
464 "use of unsupported track-status.");
465 }
466
467 if(GetProcessName() != "hElastic" && GetProcessName() != "HadronElastic"
468 && theTotalResult->GetTrackStatus()==fAlive
[962]469 && aR->GetStatusChange()==isAlive)
470 {
[819]471 // Use for debugging: G4double newWeight = theTotalResult->GetParentWeight();
472
473 G4double newKE = std::max(DBL_MIN, aR->GetEnergyChange());
[1340]474 G4DynamicParticle* aNew = new G4DynamicParticle(aT.GetParticleDefinition(),
[819]475 aR->GetMomentumChange(),
476 newKE);
477 G4HadSecondary* theSec = new G4HadSecondary(aNew, 1.0);
478 aR->AddSecondary(theSec);
479 aR->SetStatusChange(stopAndKill);
480
481 theTotalResult->ProposeTrackStatus(fStopAndKill);
482 theTotalResult->ProposeEnergy( 0.0 );
483
484 }
485 theTotalResult->ProposeLocalEnergyDeposit(aR->GetLocalEnergyDeposit());
486 theTotalResult->SetNumberOfSecondaries(aR->GetNumberOfSecondaries());
487
488 if(aR->GetStatusChange() != stopAndKill)
489 {
[1340]490 G4double newM=aT.GetParticleDefinition()->GetPDGMass();
[819]491 G4double newE=aR->GetEnergyChange() + newM;
492 G4double newP=std::sqrt(newE*newE - newM*newM);
493 G4ThreeVector newPV = newP*aR->GetMomentumChange();
494 G4LorentzVector newP4(newE, newPV);
495 newP4.rotate(rotation, it);
496 newP4*=aR->GetTrafoToLab();
497 theTotalResult->ProposeMomentumDirection(newP4.vect().unit());
498 }
499
[1340]500 for(G4int i=0; i<aR->GetNumberOfSecondaries(); ++i)
[819]501 {
502 G4LorentzVector theM = aR->GetSecondary(i)->GetParticle()->Get4Momentum();
503 theM.rotate(rotation, it);
504 theM*=aR->GetTrafoToLab();
505 aR->GetSecondary(i)->GetParticle()->Set4Momentum(theM);
506 G4double time = aR->GetSecondary(i)->GetTime();
507 if(time<0) time = aT.GetGlobalTime();
508
509 G4Track* track = new G4Track(aR->GetSecondary(i)->GetParticle(),
[962]510 time,
511 aT.GetPosition());
[819]512
513 G4double newWeight = aT.GetWeight()*aR->GetSecondary(i)->GetWeight();
514 //static G4double pinelcount=0;
[1340]515 if(xBiasOn) { newWeight *= XBiasSecondaryWeight(); }
[819]516 // G4cout << "#### ParticleDebug "
517 // <<GetProcessName()<<" "
518 // <<aR->GetSecondary(i)->GetParticle()->GetDefinition()->GetParticleName()<<" "
519 // <<aScaleFactor<<" "
520 // <<XBiasSurvivalProbability()<<" "
521 // <<XBiasSecondaryWeight()<<" "
522 // <<aT.GetWeight()<<" "
523 // <<aR->GetSecondary(i)->GetWeight()<<" "
524 // <<aR->GetSecondary(i)->GetParticle()->Get4Momentum()<<" "
525 // <<G4endl;
526 track->SetWeight(newWeight);
[962]527 track->SetTouchableHandle(aT.GetTouchableHandle());
528 theTotalResult->AddSecondary(track);
[819]529 }
530
531 aR->Clear();
532 return;
533}
[962]534
535G4IsoParticleChange* G4HadronicProcess::GetIsotopeProductionInfo()
536{
537 G4IsoParticleChange * anIsoResult = theIsoResult;
538 if(theIsoResult) theOldIsoResult = theIsoResult;
539 theIsoResult = 0;
540 return anIsoResult;
541}
542
543void G4HadronicProcess::BiasCrossSectionByFactor(G4double aScale)
544{
545 xBiasOn = true;
546 aScaleFactor = aScale;
547 G4String it = GetProcessName();
548 if( (it != "PhotonInelastic") &&
549 (it != "ElectroNuclear") &&
550 (it != "PositronNuclear") )
551 {
[1315]552 G4Exception("G4HadronicProcess::BiasCrossSectionByFactor", "007", FatalException,
[962]553 "Cross-section biasing available only for gamma and electro nuclear reactions.");
554 }
555 if(aScale<100)
556 {
[1315]557 G4Exception("G4HadronicProcess::BiasCrossSectionByFactor", "001", JustWarning,
[962]558 "Cross-section bias readjusted to be above safe limit. New value is 100");
559 aScaleFactor = 100.;
560 }
561}
[1315]562
563void
564G4HadronicProcess::CheckEnergyMomentumConservation(const G4Track& aTrack,
565 const G4Nucleus& aNucleus)
566{
567 G4double targetMass = G4NucleiProperties::GetNuclearMass(aNucleus.GetN(), aNucleus.GetZ());
568 G4LorentzVector projectile4mom = aTrack.GetDynamicParticle()->Get4Momentum();
569 G4LorentzVector target4mom(0, 0, 0, targetMass);
570 G4LorentzVector initial4mom = projectile4mom + target4mom;
571
572 G4Track* sec;
573 G4LorentzVector final4mom;
574 G4int nSec = theTotalResult->GetNumberOfSecondaries();
575 for (G4int i = 0; i < nSec; i++) {
576 sec = theTotalResult->GetSecondary(i);
577 final4mom += sec->GetDynamicParticle()->Get4Momentum();
578 }
579
580 G4LorentzVector diff = initial4mom - final4mom;
581 G4double absolute = diff.e();
582 G4double relative = absolute/aTrack.GetKineticEnergy();
583
584 G4String processName = GetProcessName();
585 G4HadronicInteraction* theModel = GetHadronicInteraction();
586 G4String modelName("none");
587 if (theModel) modelName = theModel->GetModelName();
588
589 std::pair<G4double, G4double> checkLevels = epCheckLevels;;
590 if (!levelsSetByProcess) {
591 if (theModel) checkLevels = theModel->GetEnergyMomentumCheckLevels();
592 }
593
594 G4bool relPass = false;
595 G4String relResult = "fail";
596 if (std::abs(relative) < checkLevels.first) {
597 relPass = true;
598 relResult = "pass";
599 }
600
601 G4bool absPass = false;
602 G4String absResult = "fail";
603 if (std::abs(absolute) < checkLevels.second) {
604 absPass = true;
605 absResult = "pass";
606 }
607
608 // Options for level of reporting detail:
609 // 0. off
610 // 1. report only when E/p not conserved
611 // 2. report regardless of E/p conservation
612 // 3. report only when E/p not conserved, with model names, process names, and limits
613 // 4. report regardless of E/p conservation, with model names, process names, and limits
614
615 if(epReportLevel == 4) {
616 G4cout << " Process: " << processName << " , Model: " << modelName << G4endl;
617 G4cout << " relative limit " << checkLevels.first << " relative value = "
618 << relative << " " << relResult << G4endl;
619 G4cout << " absolute limit " << checkLevels.second << " absolute value = "
620 << absolute << " " << absResult << G4endl;
621
622 } else if(epReportLevel == 3) {
623 if (!absPass || !relPass) {
624 G4cout << " Process: " << processName << " , Model: " << modelName << G4endl;
625 G4cout << " relative limit " << checkLevels.first << " relative value = "
626 << relative << " " << relResult << G4endl;
627 G4cout << " absolute limit " << checkLevels.second << " absolute value = "
628 << absolute << " " << absResult << G4endl;
629 }
630
631 } else if(epReportLevel == 2) {
632 G4cout << " relative value = " << relative << " " << relPass
633 << " absolute value = " << absolute << " " << absPass << G4endl;
634
635 } else if(epReportLevel == 1) {
636 if (!absPass || !relPass) {
637 G4cout << " relative value = " << relative << " " << relPass
638 << " absolute value = " << absolute << " " << absPass << G4endl;
639 }
640 }
641}
642
643
644void G4HadronicProcess::DumpState(const G4Track& aTrack, const G4String& method)
645{
646 G4cout << "Unrecoverable error in the method " << method << " of "
647 << GetProcessName() << G4endl;
648 G4cout << "TrackID= "<< aTrack.GetTrackID() << " ParentID= " << aTrack.GetParentID()
[1340]649 << " " << aTrack.GetParticleDefinition()->GetParticleName() << G4endl;
[1315]650 G4cout << "Ekin(GeV)= " << aTrack.GetKineticEnergy()/CLHEP::GeV
651 << "; direction= " << aTrack.GetMomentumDirection() << G4endl;
652 G4cout << "Position(mm)= " << aTrack.GetPosition()/CLHEP::mm
653 << "; material " << aTrack.GetMaterial()->GetName() << G4endl;
654 G4cout << "PhysicalVolume <" << aTrack.GetVolume()->GetName() << ">" << G4endl;
655}
656
[819]657/* end of file */
Note: See TracBrowser for help on using the repository browser.