Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/biasing/GB04/src/GB04BOptnBremSplitting.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 //
 27 /// \file GB04BOptnBremSplitting.cc
 28 /// \brief Implementation of the GB04BOptnBremSplitting class
 29 
 30 #include "GB04BOptnBremSplitting.hh"
 31 
 32 #include "G4BiasingProcessInterface.hh"
 33 #include "G4ParticleChangeForLoss.hh"
 34 
 35 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 36 
 37 GB04BOptnBremSplitting::GB04BOptnBremSplitting(G4String name)
 38   : G4VBiasingOperation(name), fSplittingFactor(1), fParticleChange()
 39 {}
 40 
 41 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 42 
 43 GB04BOptnBremSplitting::~GB04BOptnBremSplitting() {}
 44 
 45 G4VParticleChange*
 46 GB04BOptnBremSplitting::ApplyFinalStateBiasing(const G4BiasingProcessInterface* callingProcess,
 47                                                const G4Track* track, const G4Step* step, G4bool&)
 48 {
 49   // -- Collect brem. process (wrapped process) final state:
 50   G4VParticleChange* processFinalState =
 51     callingProcess->GetWrappedProcess()->PostStepDoIt(*track, *step);
 52 
 53   // -- if no splitting requested, let the brem. process to return directly its
 54   // -- generated final state:
 55   if (fSplittingFactor == 1) return processFinalState;
 56 
 57   // -- a special case here: the brem. process corrects for cross-section change
 58   // -- over the step due to energy loss by sometimes "abandoning" the interaction,
 59   // -- returning an unchanged incoming electron/positron.
 60   // -- We respect this correction, and if no secondary is produced, its means this
 61   // -- case is happening:
 62   if (processFinalState->GetNumberOfSecondaries() == 0) return processFinalState;
 63 
 64   // -- Now start the biasing:
 65   // --   - the electron state will be taken as the first one produced by the brem.
 66   // --     process, hence the one stored in above processFinalState particle change.
 67   // --     This state will be stored in our fParticleChange object.
 68   // --   - the photon accompagnying the electron will be stored also this way.
 69   // --   - we will then do fSplittingFactor - 1 call to the brem. process to collect
 70   // --     fSplittingFactor - 1 additionnal gammas. All these will be stored in our
 71   // --     fParticleChange object.
 72 
 73   // -- We called the brem. process above. Its concrete particle change is indeed
 74   // -- a "G4ParticleChangeForLoss" object. We cast this particle change to access
 75   // -- methods of the concrete G4ParticleChangeForLoss type:
 76   G4ParticleChangeForLoss* actualParticleChange = (G4ParticleChangeForLoss*)processFinalState;
 77 
 78   fParticleChange.Initialize(*track);
 79 
 80   // -- Store electron final state:
 81   fParticleChange.ProposeTrackStatus(actualParticleChange->GetTrackStatus());
 82   fParticleChange.ProposeEnergy(actualParticleChange->GetProposedKineticEnergy());
 83   fParticleChange.ProposeMomentumDirection(actualParticleChange->GetProposedMomentumDirection());
 84 
 85   // -- Now deal with the gamma's:
 86   // -- their common weight:
 87   G4double gammaWeight = track->GetWeight() / fSplittingFactor;
 88 
 89   // -- inform we will have fSplittingFactor gamma's:
 90   fParticleChange.SetNumberOfSecondaries(fSplittingFactor);
 91 
 92   // -- inform we take care of secondaries weight (otherwise these
 93   // -- secondaries are by default given the primary weight).
 94   fParticleChange.SetSecondaryWeightByProcess(true);
 95 
 96   // -- Store first gamma:
 97   G4Track* gammaTrack = actualParticleChange->GetSecondary(0);
 98   gammaTrack->SetWeight(gammaWeight);
 99   fParticleChange.AddSecondary(gammaTrack);
100   // -- and clean-up the brem. process particle change:
101   actualParticleChange->Clear();
102 
103   // -- now start the fSplittingFactor-1 calls to the brem. process to store each
104   // -- related gamma:
105   G4int nCalls = 1;
106   while (nCalls < fSplittingFactor) {
107     // ( note: we don't need to cast to actual type here, as methods for accessing
108     //   secondary particles are from base class G4VParticleChange )
109     processFinalState = callingProcess->GetWrappedProcess()->PostStepDoIt(*track, *step);
110     if (processFinalState->GetNumberOfSecondaries() == 1) {
111       gammaTrack = processFinalState->GetSecondary(0);
112       gammaTrack->SetWeight(gammaWeight);
113       fParticleChange.AddSecondary(gammaTrack);
114       nCalls++;
115     }
116     // -- very rare special case: we ignore for now.
117     else if (processFinalState->GetNumberOfSecondaries() > 1) {
118       for (G4int i = 0; i < processFinalState->GetNumberOfSecondaries(); i++)
119         delete processFinalState->GetSecondary(i);
120     }
121     processFinalState->Clear();
122   }
123 
124   // -- we are done:
125   return &fParticleChange;
126 }
127 
128 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
129