|
||||
File indexing completed on 2025-01-18 09:57:51
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 // Class: G4AdjointAlongStepWeightCorrection 0028 // Author: L. Desorgher 0029 // Organisation: SpaceIT GmbH 0030 ///////////////////////////////////////////////////////////////////////////////// 0031 // 0032 // Documentation: 0033 // 0034 // Continuous processes act on adjoint particles to continuously correct their 0035 // weight during the adjoint reverse tracking. This process is needed when 0036 // the adjoint cross sections are not scaled such that the total adjoint cross 0037 // section matches the total forward cross section. By default the mode where 0038 // the total adjoint cross section is equal to the total forward cross section 0039 // is used and therefore this along step weightcorrection factor is 1. However 0040 // in some cases (some energy ranges) the total forward cross section or the 0041 // total adjoint cross section can be zero. In this case the along step weight 0042 // correction is needed and is given by exp(-(Sigma_tot_adj-Sigma_tot_fwd).dx) 0043 // 0044 //------------------------------------------------------------- 0045 0046 #ifndef G4AdjointAlongStepWeightCorrection_h 0047 #define G4AdjointAlongStepWeightCorrection_h 1 0048 0049 #include "globals.hh" 0050 #include "G4VContinuousProcess.hh" 0051 0052 class G4AdjointCSManager; 0053 class G4MaterialCutsCouple; 0054 class G4ParticleDefinition; 0055 class G4ParticleChange; 0056 class G4Step; 0057 class G4Track; 0058 0059 class G4AdjointAlongStepWeightCorrection : public G4VContinuousProcess 0060 { 0061 public: 0062 explicit G4AdjointAlongStepWeightCorrection( 0063 const G4String& name = "ContinuousWeightCorrection", 0064 G4ProcessType type = fElectromagnetic); 0065 0066 ~G4AdjointAlongStepWeightCorrection() override; 0067 0068 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&) override; 0069 0070 void ProcessDescription(std::ostream&) const override; 0071 void DumpInfo() const override { ProcessDescription(G4cout); }; 0072 0073 G4AdjointAlongStepWeightCorrection(G4AdjointAlongStepWeightCorrection&) = 0074 delete; 0075 G4AdjointAlongStepWeightCorrection& operator=( 0076 const G4AdjointAlongStepWeightCorrection& right) = delete; 0077 0078 protected: 0079 G4double GetContinuousStepLimit(const G4Track& track, 0080 G4double previousStepSize, 0081 G4double currentMinimumStep, 0082 G4double& currentSafety) override; 0083 0084 private: 0085 void DefineMaterial(const G4MaterialCutsCouple* couple); 0086 0087 const G4MaterialCutsCouple* fCurrentCouple = nullptr; 0088 G4AdjointCSManager* fCSManager = nullptr; 0089 G4ParticleChange* fParticleChange; 0090 0091 G4double fPreStepKinEnergy = 1.; 0092 }; 0093 0094 inline void G4AdjointAlongStepWeightCorrection::DefineMaterial( 0095 const G4MaterialCutsCouple* couple) 0096 { 0097 if(couple != fCurrentCouple) 0098 { 0099 fCurrentCouple = couple; 0100 } 0101 } 0102 0103 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |