Back to home page

EIC code displayed by LXR

 
 

    


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