Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /geant4/examples/extended/biasing/GB01/src/GB01BOptrMultiParticleChangeCrossSection.cc was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

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 GB01/src/GB01BOptrMultiParticleChangeCrossSection.cc
0027 /// \brief Implementation of the GB01BOptrMultiParticleChangeCrossSection class
0028 //
0029 #include "GB01BOptrMultiParticleChangeCrossSection.hh"
0030 
0031 #include "GB01BOptrChangeCrossSection.hh"
0032 
0033 #include "G4BiasingProcessInterface.hh"
0034 #include "G4ParticleDefinition.hh"
0035 #include "G4ParticleTable.hh"
0036 #include "G4SystemOfUnits.hh"
0037 
0038 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0039 
0040 GB01BOptrMultiParticleChangeCrossSection::GB01BOptrMultiParticleChangeCrossSection()
0041   : G4VBiasingOperator("TestManyExponentialTransform")
0042 {}
0043 
0044 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0045 
0046 void GB01BOptrMultiParticleChangeCrossSection::AddParticle(G4String particleName)
0047 {
0048   const G4ParticleDefinition* particle =
0049     G4ParticleTable::GetParticleTable()->FindParticle(particleName);
0050 
0051   if (particle == nullptr) {
0052     G4ExceptionDescription ed;
0053     ed << "Particle `" << particleName << "' not found !" << G4endl;
0054     G4Exception("GB01BOptrMultiParticleChangeCrossSection::AddParticle(...)", "exGB01.02",
0055                 JustWarning, ed);
0056     return;
0057   }
0058 
0059   auto optr = new GB01BOptrChangeCrossSection(particleName);
0060   fParticlesToBias.push_back(particle);
0061   fBOptrForParticle[particle] = optr;
0062 }
0063 
0064 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0065 
0066 G4VBiasingOperation* GB01BOptrMultiParticleChangeCrossSection::ProposeOccurenceBiasingOperation(
0067   const G4Track* track, const G4BiasingProcessInterface* callingProcess)
0068 {
0069   // -- examples of limitations imposed to apply the biasing:
0070   // -- limit application of biasing to primary particles only:
0071   if (track->GetParentID() != 0) return nullptr;
0072   // -- limit to at most 5 biased interactions:
0073   if (fnInteractions > 4) return nullptr;
0074   // -- and limit to a weight of at least 0.05:
0075   if (track->GetWeight() < 0.05) return nullptr;
0076 
0077   if (fCurrentOperator)
0078     return fCurrentOperator->GetProposedOccurenceBiasingOperation(track, callingProcess);
0079   else
0080     return nullptr;
0081 }
0082 
0083 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0084 
0085 void GB01BOptrMultiParticleChangeCrossSection::StartTracking(const G4Track* track)
0086 {
0087   // -- fetch the underneath biasing operator, if any, for the current particle type:
0088   const G4ParticleDefinition* definition = track->GetParticleDefinition();
0089   auto it = fBOptrForParticle.find(definition);
0090   fCurrentOperator = nullptr;
0091   if (it != fBOptrForParticle.end()) fCurrentOperator = (*it).second;
0092 
0093   // -- reset count for number of biased interactions:
0094   fnInteractions = 0;
0095 }
0096 
0097 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0098 
0099 void GB01BOptrMultiParticleChangeCrossSection::OperationApplied(
0100   const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
0101   G4VBiasingOperation* occurenceOperationApplied, G4double weightForOccurenceInteraction,
0102   G4VBiasingOperation* finalStateOperationApplied, const G4VParticleChange* particleChangeProduced)
0103 {
0104   // -- count number of biased interactions:
0105   fnInteractions++;
0106 
0107   // -- inform the underneath biasing operator that a biased interaction occured:
0108   if (fCurrentOperator)
0109     fCurrentOperator->ReportOperationApplied(callingProcess, biasingCase, occurenceOperationApplied,
0110                                              weightForOccurenceInteraction,
0111                                              finalStateOperationApplied, particleChangeProduced);
0112 }
0113 
0114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......