Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:12

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 // -------------------------------------------------------------------
0028 //
0029 // GEANT4 Class header file
0030 //
0031 //
0032 // File name:     G4EmBiasingManager
0033 //
0034 // Author:        Vladimir Ivanchenko
0035 //
0036 // Creation date: 28.07.2011
0037 //
0038 // Modifications:
0039 //
0040 // Class Description:
0041 //
0042 // It is a class providing step limit for forced process biasing
0043 
0044 // -------------------------------------------------------------------
0045 //
0046 
0047 #ifndef G4EmBiasingManager_h
0048 #define G4EmBiasingManager_h 1
0049 
0050 #include "globals.hh"
0051 #include "G4ParticleDefinition.hh"
0052 #include "G4DynamicParticle.hh"
0053 #include "Randomize.hh"
0054 #include <vector>
0055 
0056 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0057 
0058 class G4Region;
0059 class G4Track;
0060 class G4VEnergyLossProcess;
0061 class G4VEmModel;
0062 class G4MaterialCutsCouple;
0063 class G4ParticleChangeForLoss;
0064 class G4ParticleChangeForGamma;
0065 
0066 class G4EmBiasingManager 
0067 {
0068 public:
0069 
0070   G4EmBiasingManager();
0071 
0072   ~G4EmBiasingManager();
0073 
0074   void Initialise(const G4ParticleDefinition& part,
0075           const G4String& procName, G4int verbose);
0076 
0077   // default parameters are possible
0078   void ActivateForcedInteraction(G4double length = 0.0, 
0079                  const G4String& r = "");
0080 
0081   // no default parameters
0082   void ActivateSecondaryBiasing(const G4String& region, G4double factor,
0083                 G4double energyLimit);
0084 
0085   // return forced step limit
0086   G4double GetStepLimit(G4int coupleIdx, G4double previousStep);
0087 
0088   // return weight of splitting or Russian roulette
0089   // G4DynamicParticle may be deleted
0090   // two functions are required because of the different ParticleChange
0091   // ApplySecondaryBiasing() are wrappers
0092 
0093   // for G4VEmProcess 
0094   G4double ApplySecondaryBiasing(std::vector<G4DynamicParticle*>&, 
0095                  const G4Track& track,
0096                  G4VEmModel* currentModel,
0097                  G4ParticleChangeForGamma* pParticleChange,
0098                  G4double& eloss,
0099                  G4int coupleIdx,
0100                  G4double tcut,
0101                  G4double safety = 0.0);
0102 
0103   // for G4VEnergyLossProcess 
0104   G4double ApplySecondaryBiasing(std::vector<G4DynamicParticle*>&,
0105                  const G4Track& track,
0106                  G4VEmModel* currentModel,
0107                  G4ParticleChangeForLoss* pParticleChange,
0108                  G4double& eloss, 
0109              G4int coupleIdx,  
0110                  G4double tcut, 
0111                  G4double safety = 0.0);
0112 
0113   // for G4VEnergyLossProcess 
0114   G4double ApplySecondaryBiasing(std::vector<G4Track*>&,
0115                  G4int coupleIdx);
0116 
0117   inline G4bool SecondaryBiasingRegion(G4int coupleIdx);
0118 
0119   inline G4bool ForcedInteractionRegion(G4int coupleIdx);
0120 
0121   inline void ResetForcedInteraction();
0122 
0123   G4bool CheckDirection(G4ThreeVector pos, G4ThreeVector momdir) const;
0124 
0125   G4bool  GetDirectionalSplitting() { return fDirectionalSplitting; }
0126   void    SetDirectionalSplitting(G4bool v) { fDirectionalSplitting = v; }
0127 
0128   void SetDirectionalSplittingTarget(G4ThreeVector v) 
0129     { fDirectionalSplittingTarget = v; }
0130   void SetDirectionalSplittingRadius(G4double r) 
0131     { fDirectionalSplittingRadius = r; }
0132   G4double GetWeight(G4int i);
0133 
0134   // hide copy constructor and assignment operator
0135   G4EmBiasingManager(G4EmBiasingManager &) = delete;
0136   G4EmBiasingManager & operator=(const G4EmBiasingManager &right) = delete;
0137 
0138 private:
0139 
0140   void ApplyRangeCut(std::vector<G4DynamicParticle*>& vd,
0141              const G4Track& track,
0142              G4double& eloss, 
0143              G4double safety);
0144 
0145   G4double ApplySplitting(std::vector<G4DynamicParticle*>& vd,
0146               const G4Track& track,
0147               G4VEmModel* currentModel, 
0148               G4int index,
0149               G4double tcut);
0150 
0151   G4double ApplyDirectionalSplitting(std::vector<G4DynamicParticle*>& vd,
0152         const G4Track& track,
0153         G4VEmModel* currentModel, 
0154         G4int index,
0155         G4double tcut,
0156         G4ParticleChangeForGamma* partChange);
0157 
0158   G4double ApplyDirectionalSplitting(std::vector<G4DynamicParticle*>& vd,
0159         const G4Track& track,
0160         G4VEmModel* currentModel, 
0161         G4int index,
0162         G4double tcut);
0163 
0164   inline G4double ApplyRussianRoulette(std::vector<G4DynamicParticle*>& vd,
0165                        G4int index); 
0166 
0167   G4VEnergyLossProcess* eIonisation = nullptr;
0168 
0169   const G4ParticleDefinition* theElectron;
0170   const G4ParticleDefinition* theGamma;
0171 
0172   G4double fSafetyMin;
0173   G4double currentStepLimit = 0.0;
0174   G4double fDirectionalSplittingRadius = 0.0;
0175 
0176   G4int nForcedRegions = 0;
0177   G4int nSecBiasedRegions = 0;
0178 
0179   G4bool startTracking = true;
0180   G4bool fDirectionalSplitting = false;
0181 
0182   G4ThreeVector fDirectionalSplittingTarget;
0183   std::vector<G4double> fDirectionalSplittingWeights;
0184 
0185   std::vector<G4double>        lengthForRegion;
0186   std::vector<G4double>        secBiasedWeight;
0187   std::vector<G4double>        secBiasedEnegryLimit;
0188 
0189   std::vector<const G4Region*> forcedRegions;
0190   std::vector<const G4Region*> secBiasedRegions;
0191 
0192   std::vector<G4int>           nBremSplitting;
0193   std::vector<G4int>           idxForcedCouple;
0194   std::vector<G4int>           idxSecBiasedCouple;
0195 
0196   std::vector<G4DynamicParticle*> tmpSecondaries;
0197 };
0198 
0199 inline G4bool 
0200 G4EmBiasingManager::SecondaryBiasingRegion(G4int coupleIdx)
0201 {
0202   G4bool res = false;
0203   if(nSecBiasedRegions > 0) {
0204     if(idxSecBiasedCouple[coupleIdx] >= 0) { res = true; }
0205   }
0206   return res;
0207 }
0208 
0209 inline G4bool G4EmBiasingManager::ForcedInteractionRegion(G4int coupleIdx)
0210 {
0211   G4bool res = false;
0212   if(nForcedRegions > 0) {
0213     if(idxForcedCouple[coupleIdx] >= 0) { res = true; }
0214   }
0215   return res;
0216 }
0217 
0218 inline void G4EmBiasingManager::ResetForcedInteraction()
0219 {
0220   startTracking = true;
0221 }
0222 
0223 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0224 
0225 inline G4double
0226 G4EmBiasingManager::ApplyRussianRoulette(std::vector<G4DynamicParticle*>& vd,
0227                      G4int index)
0228 {
0229   size_t n = vd.size();
0230   G4double weight = secBiasedWeight[index];
0231   for(size_t k=0; k<n; ++k) {
0232     if(G4UniformRand()*weight > 1.0) {
0233       const G4DynamicParticle* dp = vd[k];
0234       delete dp;
0235       vd[k] = nullptr;
0236     }
0237   }
0238   return weight;
0239 }
0240 
0241 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0242 
0243 #endif