Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-12-16 09:30:09

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 /// \file GB02/src/GB02BOptrMultiParticleForceCollision.cc
0027 /// \brief Implementation of the GB02BOptrMultiParticleForceCollision class
0028 //
0029 #include "GB02BOptrMultiParticleForceCollision.hh"
0030 
0031 #include "G4BOptrForceCollision.hh"
0032 #include "G4BiasingProcessInterface.hh"
0033 #include "G4ParticleDefinition.hh"
0034 #include "G4ParticleTable.hh"
0035 
0036 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0037 
0038 GB02BOptrMultiParticleForceCollision::GB02BOptrMultiParticleForceCollision()
0039   : G4VBiasingOperator("TestManyForceCollision")
0040 {}
0041 
0042 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0043 
0044 void GB02BOptrMultiParticleForceCollision::AddParticle(G4String particleName)
0045 {
0046   const G4ParticleDefinition* particle =
0047     G4ParticleTable::GetParticleTable()->FindParticle(particleName);
0048 
0049   if (particle == nullptr) {
0050     G4ExceptionDescription ed;
0051     ed << "Particle `" << particleName << "' not found !" << G4endl;
0052     G4Exception("GB02BOptrMultiParticleForceCollision::AddParticle(...)", "exGB02.01", JustWarning,
0053                 ed);
0054     return;
0055   }
0056 
0057   auto optr = new G4BOptrForceCollision(particleName, "ForceCollisionFor" + particleName);
0058   fParticlesToBias.push_back(particle);
0059   fBOptrForParticle[particle] = optr;
0060 }
0061 
0062 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0063 
0064 G4VBiasingOperation* GB02BOptrMultiParticleForceCollision::ProposeOccurenceBiasingOperation(
0065   const G4Track* track, const G4BiasingProcessInterface* callingProcess)
0066 {
0067   if (fCurrentOperator)
0068     return fCurrentOperator->GetProposedOccurenceBiasingOperation(track, callingProcess);
0069   else
0070     return nullptr;
0071 }
0072 
0073 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0074 
0075 G4VBiasingOperation* GB02BOptrMultiParticleForceCollision::ProposeNonPhysicsBiasingOperation(
0076   const G4Track* track, const G4BiasingProcessInterface* callingProcess)
0077 {
0078   if (fCurrentOperator)
0079     return fCurrentOperator->GetProposedNonPhysicsBiasingOperation(track, callingProcess);
0080   else
0081     return nullptr;
0082 }
0083 
0084 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0085 
0086 G4VBiasingOperation* GB02BOptrMultiParticleForceCollision::ProposeFinalStateBiasingOperation(
0087   const G4Track* track, const G4BiasingProcessInterface* callingProcess)
0088 {
0089   if (fCurrentOperator)
0090     return fCurrentOperator->GetProposedFinalStateBiasingOperation(track, callingProcess);
0091   else
0092     return nullptr;
0093 }
0094 
0095 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0096 
0097 void GB02BOptrMultiParticleForceCollision::StartTracking(const G4Track* track)
0098 {
0099   const G4ParticleDefinition* definition = track->GetParticleDefinition();
0100   auto it = fBOptrForParticle.find(definition);
0101   fCurrentOperator = nullptr;
0102   if (it != fBOptrForParticle.end()) fCurrentOperator = (*it).second;
0103 }
0104 
0105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0106 
0107 void GB02BOptrMultiParticleForceCollision::OperationApplied(
0108   const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
0109   G4VBiasingOperation* operationApplied, const G4VParticleChange* particleChangeProduced)
0110 {
0111   if (fCurrentOperator)
0112     fCurrentOperator->ReportOperationApplied(callingProcess, biasingCase, operationApplied,
0113                                              particleChangeProduced);
0114 }
0115 
0116 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0117 
0118 void GB02BOptrMultiParticleForceCollision::OperationApplied(
0119   const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
0120   G4VBiasingOperation* occurenceOperationApplied, G4double weightForOccurenceInteraction,
0121   G4VBiasingOperation* finalStateOperationApplied, const G4VParticleChange* particleChangeProduced)
0122 {
0123   if (fCurrentOperator)
0124     fCurrentOperator->ReportOperationApplied(callingProcess, biasingCase, occurenceOperationApplied,
0125                                              weightForOccurenceInteraction,
0126                                              finalStateOperationApplied, particleChangeProduced);
0127 }
0128 
0129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0130 
0131 void GB02BOptrMultiParticleForceCollision::ExitBiasing(
0132   const G4Track* track, const G4BiasingProcessInterface* callingProcess)
0133 {
0134   if (fCurrentOperator) fCurrentOperator->ExitingBiasing(track, callingProcess);
0135 }
0136 
0137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......