Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/biasing/GB02/src/GB02BOptrMultiParticleForceCollision.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 26 /// \file GB02/src/GB02BOptrMultiParticleForceCollision.cc
 27 /// \brief Implementation of the GB02BOptrMultiParticleForceCollision class
 28 //
 29 #include "GB02BOptrMultiParticleForceCollision.hh"
 30 
 31 #include "G4BOptrForceCollision.hh"
 32 #include "G4BiasingProcessInterface.hh"
 33 #include "G4ParticleDefinition.hh"
 34 #include "G4ParticleTable.hh"
 35 
 36 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 37 
 38 GB02BOptrMultiParticleForceCollision::GB02BOptrMultiParticleForceCollision()
 39   : G4VBiasingOperator("TestManyForceCollision")
 40 {}
 41 
 42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 43 
 44 void GB02BOptrMultiParticleForceCollision::AddParticle(G4String particleName)
 45 {
 46   const G4ParticleDefinition* particle =
 47     G4ParticleTable::GetParticleTable()->FindParticle(particleName);
 48 
 49   if (particle == 0) {
 50     G4ExceptionDescription ed;
 51     ed << "Particle `" << particleName << "' not found !" << G4endl;
 52     G4Exception("GB02BOptrMultiParticleForceCollision::AddParticle(...)", "exGB02.01", JustWarning,
 53                 ed);
 54     return;
 55   }
 56 
 57   G4BOptrForceCollision* optr =
 58     new G4BOptrForceCollision(particleName, "ForceCollisionFor" + particleName);
 59   fParticlesToBias.push_back(particle);
 60   fBOptrForParticle[particle] = optr;
 61 }
 62 
 63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 64 
 65 G4VBiasingOperation* GB02BOptrMultiParticleForceCollision::ProposeOccurenceBiasingOperation(
 66   const G4Track* track, const G4BiasingProcessInterface* callingProcess)
 67 {
 68   if (fCurrentOperator)
 69     return fCurrentOperator->GetProposedOccurenceBiasingOperation(track, callingProcess);
 70   else
 71     return 0;
 72 }
 73 
 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 75 
 76 G4VBiasingOperation* GB02BOptrMultiParticleForceCollision::ProposeNonPhysicsBiasingOperation(
 77   const G4Track* track, const G4BiasingProcessInterface* callingProcess)
 78 {
 79   if (fCurrentOperator)
 80     return fCurrentOperator->GetProposedNonPhysicsBiasingOperation(track, callingProcess);
 81   else
 82     return 0;
 83 }
 84 
 85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 86 
 87 G4VBiasingOperation* GB02BOptrMultiParticleForceCollision::ProposeFinalStateBiasingOperation(
 88   const G4Track* track, const G4BiasingProcessInterface* callingProcess)
 89 {
 90   if (fCurrentOperator)
 91     return fCurrentOperator->GetProposedFinalStateBiasingOperation(track, callingProcess);
 92   else
 93     return 0;
 94 }
 95 
 96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 97 
 98 void GB02BOptrMultiParticleForceCollision::StartTracking(const G4Track* track)
 99 {
100   const G4ParticleDefinition* definition = track->GetParticleDefinition();
101   std::map<const G4ParticleDefinition*, G4BOptrForceCollision*>::iterator it =
102     fBOptrForParticle.find(definition);
103   fCurrentOperator = 0;
104   if (it != fBOptrForParticle.end()) fCurrentOperator = (*it).second;
105 }
106 
107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108 
109 void GB02BOptrMultiParticleForceCollision::OperationApplied(
110   const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
111   G4VBiasingOperation* operationApplied, const G4VParticleChange* particleChangeProduced)
112 {
113   if (fCurrentOperator)
114     fCurrentOperator->ReportOperationApplied(callingProcess, biasingCase, operationApplied,
115                                              particleChangeProduced);
116 }
117 
118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
119 
120 void GB02BOptrMultiParticleForceCollision::OperationApplied(
121   const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
122   G4VBiasingOperation* occurenceOperationApplied, G4double weightForOccurenceInteraction,
123   G4VBiasingOperation* finalStateOperationApplied, const G4VParticleChange* particleChangeProduced)
124 {
125   if (fCurrentOperator)
126     fCurrentOperator->ReportOperationApplied(callingProcess, biasingCase, occurenceOperationApplied,
127                                              weightForOccurenceInteraction,
128                                              finalStateOperationApplied, particleChangeProduced);
129 }
130 
131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132 
133 void GB02BOptrMultiParticleForceCollision::ExitBiasing(
134   const G4Track* track, const G4BiasingProcessInterface* callingProcess)
135 {
136   if (fCurrentOperator) fCurrentOperator->ExitingBiasing(track, callingProcess);
137 }
138 
139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140