Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-12-18 09:32:26

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 GB05BOptnSplitAndKillByCrossSection.cc
0028 /// \brief Implementation of the GB05BOptnSplitAndKillByCrossSection class
0029 
0030 #include "GB05BOptnSplitAndKillByCrossSection.hh"
0031 
0032 #include "Randomize.hh"
0033 
0034 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0035 
0036 GB05BOptnSplitAndKillByCrossSection::GB05BOptnSplitAndKillByCrossSection(G4String name)
0037   : G4VBiasingOperation(name), fParticleChange(), fInteractionLength(-1.0)
0038 {}
0039 
0040 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0041 
0042 GB05BOptnSplitAndKillByCrossSection::~GB05BOptnSplitAndKillByCrossSection() = default;
0043 
0044 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0045 
0046 G4double GB05BOptnSplitAndKillByCrossSection::DistanceToApplyOperation(const G4Track*, G4double,
0047                                                                        G4ForceCondition* condition)
0048 {
0049   *condition = NotForced;
0050 
0051   // -- Sample the exponential law using the total interaction length of processes
0052   // -- to couterbalance for:
0053   G4double proposedStepLength = -std::log(G4UniformRand()) * fInteractionLength;
0054 
0055   return proposedStepLength;
0056 }
0057 
0058 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0059 
0060 G4VParticleChange*
0061 GB05BOptnSplitAndKillByCrossSection::GenerateBiasingFinalState(const G4Track* track, const G4Step*)
0062 {
0063   // -- This method is called if we have limited the step.
0064   // -- We hence make the splitting or killing.
0065 
0066   // Get track weight:
0067   G4double initialWeight = track->GetWeight();
0068 
0069   // The "particle change" is the object to be used to communicate to
0070   // the tracking the update of the primary state and/or creation
0071   // secondary tracks.
0072   fParticleChange.Initialize(*track);
0073 
0074   // -- Splitting and killing factors.
0075   // -- They are taken the same, but the killing factor can be make bigger.
0076   G4double splittingFactor = 2.0;
0077   G4double killingFactor = 2.0;
0078 
0079   if (track->GetMomentumDirection().z() > 0) {
0080     // -- We split if the track is moving forward:
0081 
0082     // Define the tracks weight:
0083     G4double weightOfTrack = initialWeight / splittingFactor;
0084 
0085     // Ask currect track weight to be changed to new value:
0086     fParticleChange.ProposeParentWeight(weightOfTrack);
0087     // Now make clones of this track (this is the actual splitting):
0088     // we will then have the primary and clone of it, hence the
0089     // splitting by a factor 2:
0090     auto clone = new G4Track(*track);
0091     clone->SetWeight(weightOfTrack);
0092     fParticleChange.AddSecondary(clone);
0093     // -- Below's call added for safety & illustration : inform particle change to not
0094     // -- modify the clone (ie : daughter) weight to male it that of the
0095     // -- primary. Here call not mandatory and both tracks have same weights.
0096     fParticleChange.SetSecondaryWeightByProcess(true);
0097   }
0098   else {
0099     // -- We apply Russian roulette if the track is moving backward:
0100 
0101     // Shoot a random number (in ]0,1[ segment):
0102     G4double random = G4UniformRand();
0103     G4double killingProbability = 1.0 - 1.0 / killingFactor;
0104     if (random < killingProbability) {
0105       // We ask for the the track to be killed:
0106       fParticleChange.ProposeTrackStatus(fStopAndKill);
0107     }
0108     else {
0109       // In this case, the track survives. We change its weight
0110       // to conserve weight among killed and survival tracks:
0111       fParticleChange.ProposeParentWeight(initialWeight * killingFactor);
0112     }
0113   }
0114 
0115   return &fParticleChange;
0116 }
0117 
0118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......