Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:20:45

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 GB04BOptnBremSplitting.cc
0028 /// \brief Implementation of the GB04BOptnBremSplitting class
0029 
0030 #include "GB04BOptnBremSplitting.hh"
0031 
0032 #include "G4BiasingProcessInterface.hh"
0033 #include "G4ParticleChangeForLoss.hh"
0034 
0035 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0036 
0037 GB04BOptnBremSplitting::GB04BOptnBremSplitting(G4String name)
0038   : G4VBiasingOperation(name), fSplittingFactor(1), fParticleChange()
0039 {}
0040 
0041 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0042 
0043 GB04BOptnBremSplitting::~GB04BOptnBremSplitting() {}
0044 
0045 G4VParticleChange*
0046 GB04BOptnBremSplitting::ApplyFinalStateBiasing(const G4BiasingProcessInterface* callingProcess,
0047                                                const G4Track* track, const G4Step* step, G4bool&)
0048 {
0049   // -- Collect brem. process (wrapped process) final state:
0050   G4VParticleChange* processFinalState =
0051     callingProcess->GetWrappedProcess()->PostStepDoIt(*track, *step);
0052 
0053   // -- if no splitting requested, let the brem. process to return directly its
0054   // -- generated final state:
0055   if (fSplittingFactor == 1) return processFinalState;
0056 
0057   // -- a special case here: the brem. process corrects for cross-section change
0058   // -- over the step due to energy loss by sometimes "abandoning" the interaction,
0059   // -- returning an unchanged incoming electron/positron.
0060   // -- We respect this correction, and if no secondary is produced, its means this
0061   // -- case is happening:
0062   if (processFinalState->GetNumberOfSecondaries() == 0) return processFinalState;
0063 
0064   // -- Now start the biasing:
0065   // --   - the electron state will be taken as the first one produced by the brem.
0066   // --     process, hence the one stored in above processFinalState particle change.
0067   // --     This state will be stored in our fParticleChange object.
0068   // --   - the photon accompagnying the electron will be stored also this way.
0069   // --   - we will then do fSplittingFactor - 1 call to the brem. process to collect
0070   // --     fSplittingFactor - 1 additionnal gammas. All these will be stored in our
0071   // --     fParticleChange object.
0072 
0073   // -- We called the brem. process above. Its concrete particle change is indeed
0074   // -- a "G4ParticleChangeForLoss" object. We cast this particle change to access
0075   // -- methods of the concrete G4ParticleChangeForLoss type:
0076   G4ParticleChangeForLoss* actualParticleChange = (G4ParticleChangeForLoss*)processFinalState;
0077 
0078   fParticleChange.Initialize(*track);
0079 
0080   // -- Store electron final state:
0081   fParticleChange.ProposeTrackStatus(actualParticleChange->GetTrackStatus());
0082   fParticleChange.ProposeEnergy(actualParticleChange->GetProposedKineticEnergy());
0083   fParticleChange.ProposeMomentumDirection(actualParticleChange->GetProposedMomentumDirection());
0084 
0085   // -- Now deal with the gamma's:
0086   // -- their common weight:
0087   G4double gammaWeight = track->GetWeight() / fSplittingFactor;
0088 
0089   // -- inform we will have fSplittingFactor gamma's:
0090   fParticleChange.SetNumberOfSecondaries(fSplittingFactor);
0091 
0092   // -- inform we take care of secondaries weight (otherwise these
0093   // -- secondaries are by default given the primary weight).
0094   fParticleChange.SetSecondaryWeightByProcess(true);
0095 
0096   // -- Store first gamma:
0097   G4Track* gammaTrack = actualParticleChange->GetSecondary(0);
0098   gammaTrack->SetWeight(gammaWeight);
0099   fParticleChange.AddSecondary(gammaTrack);
0100   // -- and clean-up the brem. process particle change:
0101   actualParticleChange->Clear();
0102 
0103   // -- now start the fSplittingFactor-1 calls to the brem. process to store each
0104   // -- related gamma:
0105   G4int nCalls = 1;
0106   while (nCalls < fSplittingFactor) {
0107     // ( note: we don't need to cast to actual type here, as methods for accessing
0108     //   secondary particles are from base class G4VParticleChange )
0109     processFinalState = callingProcess->GetWrappedProcess()->PostStepDoIt(*track, *step);
0110     if (processFinalState->GetNumberOfSecondaries() == 1) {
0111       gammaTrack = processFinalState->GetSecondary(0);
0112       gammaTrack->SetWeight(gammaWeight);
0113       fParticleChange.AddSecondary(gammaTrack);
0114       nCalls++;
0115     }
0116     // -- very rare special case: we ignore for now.
0117     else if (processFinalState->GetNumberOfSecondaries() > 1) {
0118       for (G4int i = 0; i < processFinalState->GetNumberOfSecondaries(); i++)
0119         delete processFinalState->GetSecondary(i);
0120     }
0121     processFinalState->Clear();
0122   }
0123 
0124   // -- we are done:
0125   return &fParticleChange;
0126 }
0127 
0128 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......