Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-03-31 07:49:46

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 GB05BOptrSplitAndKillByCrossSection.cc
0027 /// \brief Implementation of the GB05BOptrSplitAndKillByCrossSection class
0028 
0029 #include "GB05BOptrSplitAndKillByCrossSection.hh"
0030 
0031 #include "GB05BOptnSplitAndKillByCrossSection.hh"
0032 
0033 #include "G4BiasingProcessInterface.hh"
0034 #include "G4ParticleDefinition.hh"
0035 #include "G4ParticleTable.hh"
0036 #include "G4VProcess.hh"
0037 
0038 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0039 
0040 GB05BOptrSplitAndKillByCrossSection::GB05BOptrSplitAndKillByCrossSection(G4String particleName,
0041                                                                          G4String name)
0042   : G4VBiasingOperator(name), fSetup(true)
0043 {
0044   fParticleToBias = G4ParticleTable::GetParticleTable()->FindParticle(particleName);
0045 
0046   if (fParticleToBias == nullptr) {
0047     G4ExceptionDescription ed;
0048     ed << "Particle `" << particleName << "' not found !" << G4endl;
0049     G4Exception("GB05BOptrSplitAndKillByCrossSection(...)", "exGB05.01", JustWarning, ed);
0050   }
0051 
0052   fSplitAndKillByCrossSection =
0053     new GB05BOptnSplitAndKillByCrossSection("splitterFor_" + particleName);
0054 }
0055 
0056 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0057 
0058 GB05BOptrSplitAndKillByCrossSection::~GB05BOptrSplitAndKillByCrossSection()
0059 {
0060   delete fSplitAndKillByCrossSection;
0061 }
0062 
0063 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0064 
0065 void GB05BOptrSplitAndKillByCrossSection::StartRun()
0066 {
0067   // ---------------
0068   // -- Setup stage:
0069   // ---------------
0070   // -- Start by collecting the pointer of the physics processes
0071   // -- considered for the splitting by cross-sections. Doing so,
0072   // -- this also verifies that these physics processes are each
0073   // -- under control of a G4BiasingProcessInterface wrapper.
0074   if (fSetup) {
0075     const G4ProcessManager* processManager = fParticleToBias->GetProcessManager();
0076     const G4BiasingProcessSharedData* sharedData =
0077       G4BiasingProcessInterface::GetSharedData(processManager);
0078     if (sharedData) {
0079       for (size_t i = 0; i < fProcessesToEquipoise.size(); i++) {
0080         G4bool processFound(false);
0081         for (size_t j = 0; j < (sharedData->GetPhysicsBiasingProcessInterfaces()).size(); j++) {
0082           const G4BiasingProcessInterface* wrapperProcess =
0083             (sharedData->GetPhysicsBiasingProcessInterfaces())[j];
0084           if (fProcessesToEquipoise[i] == wrapperProcess->GetWrappedProcess()->GetProcessName()) {
0085             fProcesses.push_back(wrapperProcess->GetWrappedProcess());
0086             processFound = true;
0087             break;
0088           }
0089         }
0090         if (!processFound) {
0091           G4String particleName = "(unknown)";
0092           if (fParticleToBias != nullptr) {
0093             particleName = fParticleToBias->GetParticleName();
0094           }
0095           G4ExceptionDescription ed;
0096           ed << "Process `" << fProcessesToEquipoise[i] << "' not found for particle `"
0097              << particleName << "'" << G4endl;
0098           G4Exception("GB05BOptrSplitAndKillByCrossSection::StartRun(...)", "exGB05.02",
0099                       JustWarning, ed);
0100         }
0101       }
0102     }
0103     fSetup = false;
0104   }
0105 
0106   if (fProcessesToEquipoise.size() == 0 || fProcesses.size() == 0) {
0107     G4ExceptionDescription ed;
0108     ed << "No processes to counterbalance for defined or found ! "
0109        << "Biasing will do nothing." << G4endl;
0110     G4Exception("GB05BOptrSplitAndKillByCrossSection::StartRun(...)", "exGB05.03", JustWarning, ed);
0111   }
0112 }
0113 
0114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0115 
0116 G4VBiasingOperation* GB05BOptrSplitAndKillByCrossSection::ProposeNonPhysicsBiasingOperation(
0117   const G4Track* track, const G4BiasingProcessInterface* /*callingProcess*/)
0118 {
0119   // -----------------------------------------------------
0120   // -- Check if current particle type is the one to bias:
0121   // -----------------------------------------------------
0122   if (track->GetDefinition() != fParticleToBias) return nullptr;
0123 
0124   // --------------------------------------------------------------------
0125   // -- Compute the total cross-section for the physics processes
0126   // -- considered.
0127   // -- These physics processes have been updated to the current track
0128   // -- state by their related wrapper G4BiasingProcessInterface objects,
0129   // -- the cross-sections/mean free pathes are hence usable.
0130   // --------------------------------------------------------------------
0131   G4double totalCrossSection(0.0);
0132   for (size_t i = 0; i < fProcesses.size(); i++) {
0133     G4double interactionLength = fProcesses[i]->GetCurrentInteractionLength();
0134     if (interactionLength < DBL_MAX / 10.) totalCrossSection += 1. / interactionLength;
0135   }
0136   if (totalCrossSection < DBL_MIN) return nullptr;
0137 
0138   G4double totalInteractionLength = 1. / totalCrossSection;
0139 
0140   // ---------------------------------------------------------------------
0141   // -- Passes the updated "absorption" cross-section (interaction length)
0142   // -- to the biasing operation, and returns this operation:
0143   // ---------------------------------------------------------------------
0144   fSplitAndKillByCrossSection->SetInteractionLength(totalInteractionLength);
0145 
0146   return fSplitAndKillByCrossSection;
0147 }
0148 
0149 void GB05BOptrSplitAndKillByCrossSection::AddProcessToEquipoise(G4String processName)
0150 {
0151   fProcessesToEquipoise.push_back(processName);
0152 }