Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4BiasingOperator.cxx
Go to the documentation of this file.
1//------------------------------------------------
2// The Geant4 Virtual Monte Carlo package
3// Copyright (C) 2007 - 2019 Ivana Hrivnacova
4// All rights reserved.
5//
6// For the licensing terms see geant4_vmc/LICENSE.
7// Contact: root-vmc@cern.ch
8//-------------------------------------------------
9
14
15#include "TG4BiasingOperator.h"
16#include "TG4BiasingOperation.h"
17
18#include "G4BiasingProcessInterface.hh"
19#include "G4ParticleDefinition.hh"
20#include "G4ParticleTable.hh"
21
23{
24 fBiasingOperation = new TG4BiasingOperation("BiasingOperation");
25}
26
27void TG4BiasingOperator::AddParticle(G4String particleName)
28{
29 const G4ParticleDefinition* particle =
30 G4ParticleTable::GetParticleTable()->FindParticle(particleName);
31 if (particle == 0) {
32 G4ExceptionDescription ed;
33 ed << "Particle `" << particleName << "' not found !" << G4endl;
34 G4Exception(
35 "TG4BiasingOperator::AddParticle(...)", "BiasError", JustWarning, ed);
36 return;
37 }
38 fParticlesToBias.push_back(particle);
39}
40
42 const G4Track*, const G4BiasingProcessInterface* callingProcess)
43{
44
45 // G4cout << "In TG4BiasingOperator::ProposeFinalStateTG4BiasingOperation "
46 // << " calling process: " << callingProcess
47 // << " wrapper process: " << callingProcess->GetWrappedProcess();
48 // if ( callingProcess->GetWrappedProcess() ) {
49 // G4cout << " wrapper process name: " <<
50 // callingProcess->GetWrappedProcess()->GetProcessName();
51 // }
52 // G4cout << G4endl;
53
54 // Apply the biasing operation only for inelastic processes of:
55 // proton, neutron, pion+ and pion-
56 if (callingProcess && callingProcess->GetWrappedProcess() &&
57 (callingProcess->GetWrappedProcess()->GetProcessName() ==
58 "protonInelastic" ||
59 callingProcess->GetWrappedProcess()->GetProcessName() ==
60 "neutronInelastic" ||
61 callingProcess->GetWrappedProcess()->GetProcessName() ==
62 "pi+Inelastic" ||
63 callingProcess->GetWrappedProcess()->GetProcessName() ==
64 "pi-Inelastic")) {
65 // G4cout << "In TG4BiasingOperator: Returning " << fBiasingOperation <<
66 // G4endl;
67 return fBiasingOperation;
68 }
69 else {
70 // G4cout << "In TG4BiasingOperator: Returning 0 " << G4endl;
71 return 0;
72 }
73}
Definition of the TG4BiasingOperation class.
Definition of the TG4BiasingOperator class.
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
TG4BiasingOperation * fBiasingOperation
std::vector< const G4ParticleDefinition * > fParticlesToBias
void AddParticle(G4String particleName)