|
||||
File indexing completed on 2025-01-18 09:58:07
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 // Author: Mathieu Karamitros 0028 0029 // The code is developed in the framework of the ESA AO7146 0030 // 0031 // We would be very happy hearing from you, send us your feedback! :) 0032 // 0033 // In order for Geant4-DNA to be maintained and still open-source, 0034 // article citations are crucial. 0035 // If you use Geant4-DNA chemistry and you publish papers about your software, 0036 // in addition to the general paper on Geant4-DNA: 0037 // 0038 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178 0039 // 0040 // we would be very happy if you could please also cite the following 0041 // reference papers on chemistry: 0042 // 0043 // J. Comput. Phys. 274 (2014) 841-882 0044 // Prog. Nucl. Sci. Tec. 2 (2011) 503-508 0045 0046 #ifndef G4ITBROWNIANTRANSPORTATION_H 0047 #define G4ITBROWNIANTRANSPORTATION_H 0048 0049 #include "G4ITTransportation.hh" 0050 0051 class G4SafetyHelper; 0052 class G4Molecule; 0053 class G4VUserBrownianAction; 0054 0055 // experimental 0056 class G4BrownianAction 0057 { 0058 public: 0059 G4BrownianAction()= default; 0060 virtual ~G4BrownianAction()= default; 0061 0062 // virtual G4double GetDiffusionCoefficient(G4Material*, 0063 // G4Molecule*) { return 0;} 0064 0065 // If returns true: track is killed 0066 virtual void Transport(const G4Track&, 0067 G4ParticleChangeForTransport&) = 0; 0068 }; 0069 0070 0071 /* \brief {The transportation method implemented is the one from 0072 * Ermak-McCammon : J. Chem. Phys. 69, 1352 (1978). 0073 * To compute time and space intervals to reach a volume boundary, 0074 * there are two alternative methods proposed by this process. 0075 * 0076 * ** Method 1 selects a minimum distance to the next 0077 * boundary using to the following formula: 0078 * 0079 * --> t_min = (safety* safety) / (8 * diffusionCoefficient); 0080 * this corresponds to 5% probability of the Brownian particle to cross 0081 * the boundary - isotropic distance to nearest boundary (safety) is used 0082 * 0083 * OR if the flag "speed me up" is on: 0084 * 0085 * --> t_min = (geometryStepLength * geometryStepLength) / diffusionCoefficient; 0086 * this corresponds to 50% probability of the Brownian particle to cross 0087 * the boundary - distance along current direction to nearest boundary is used 0088 * 0089 * NB: By default, method 1 with the flag "speed me up is used". 0090 * In addition, one may want to used the minimum time step limit defined 0091 * in G4Scheduler through the G4UserTimeStepAction. If so, speed level might 0092 * be set to 2. But minimum time steps have to be set in the user class. 0093 * 0094 * ** Method 2 can randomly compute the time to the next boundary using the 0095 * following formula: 0096 * 0097 * t_random = 1 / (4 * diffusionCoefficient)* std::pow(geometryStepLength / 0098 * InvErfc(G4UniformRand()),2); 0099 * For release 10.1, this is using the 1D cumulative density function. 0100 * 0101 * At each diffusion step, the direction of the particle is randomly selected. 0102 * For now, the geometryStepLength corresponds to the distance to the 0103 * nearest boundary along the direction of diffusion which selected randomly. 0104 * 0105 * Method 2 is currently deactivated by default. 0106 * } 0107 */ 0108 0109 class G4DNABrownianTransportation : public G4ITTransportation 0110 { 0111 public: 0112 G4DNABrownianTransportation(const G4String& aName = 0113 "DNABrownianTransportation", 0114 G4int verbosityLevel = 0); 0115 ~G4DNABrownianTransportation() override; 0116 G4DNABrownianTransportation(const G4DNABrownianTransportation&) = delete; 0117 G4DNABrownianTransportation& operator=(const G4DNABrownianTransportation&) = delete; 0118 0119 inline void SetBrownianAction(G4BrownianAction*); 0120 inline void SetUserBrownianAction(G4VUserBrownianAction*); 0121 0122 void BuildPhysicsTable(const G4ParticleDefinition&) override; 0123 0124 void StartTracking(G4Track* aTrack) override; 0125 0126 void ComputeStep(const G4Track&, 0127 const G4Step&, 0128 const G4double, 0129 G4double&) override; 0130 0131 G4double 0132 AlongStepGetPhysicalInteractionLength(const G4Track& /*track*/, 0133 G4double /*previousStepSize*/, 0134 G4double /*currentMinimumStep*/, 0135 G4double& /*currentSafety*/, 0136 G4GPILSelection* /*selection*/) override; 0137 0138 G4VParticleChange* PostStepDoIt(const G4Track& track, const G4Step&) override; 0139 0140 G4VParticleChange* AlongStepDoIt(const G4Track& track, const G4Step&) override; 0141 0142 // Boundary is crossed at time at which: 0143 // * either 5% of the distribution might be over boundary - the last position 0144 // is adjusted on boundary 0145 // * or if speedUp (from level 1) is activated - 50% of the distribution might 0146 // be over boundary, the particles are also allowed to jump over boundary 0147 inline void UseMaximumTimeBeforeReachingBoundary(bool flag = true) 0148 { 0149 fUseMaximumTimeBeforeReachingBoundary = flag; 0150 } 0151 0152 // Random sampling time at which boundary is crossed 0153 // WARNING: For release 10.1, this is a 1D approximation for sampling time 0154 // but 3D for diffusion 0155 // If speed up IS activated, particles are allowed jump over barrier 0156 inline void UseCumulativeDensitFunction(bool flag = true) 0157 { 0158 if(flag) 0159 { 0160 fUseMaximumTimeBeforeReachingBoundary = false; 0161 return; 0162 } 0163 fUseMaximumTimeBeforeReachingBoundary = true; 0164 } 0165 0166 // Use limiting time steps defined in the scheduler 0167 inline void UseLimitingTimeSteps(bool flag = true) 0168 { 0169 fUseSchedulerMinTimeSteps = flag; 0170 } 0171 0172 inline void SpeedLevel(int level) 0173 { 0174 if(level < 0) level =0; 0175 else if(level > 2) level = 2; 0176 0177 switch(level) 0178 { 0179 case 0: 0180 fSpeedMeUp = false; 0181 fUseSchedulerMinTimeSteps = false; 0182 return; 0183 0184 case 1: 0185 fSpeedMeUp = true; 0186 fUseSchedulerMinTimeSteps = false; 0187 return; 0188 0189 case 2: 0190 //====================================================================== 0191 // NB: BE AWARE THAT IF NO MIN TIME STEPS NO TIME STEPS HAVE BEEN 0192 // PROVIDED TO G4Scheduler THIS LEVEL MIGHT BE SLOWER THAN LEVEL 1 0193 //====================================================================== 0194 fSpeedMeUp = true; 0195 fUseSchedulerMinTimeSteps = true; 0196 return; 0197 } 0198 } 0199 0200 protected: 0201 0202 G4double ComputeGeomLimit(const G4Track& track, 0203 G4double& presafety, 0204 G4double limit); 0205 0206 void Diffusion(const G4Track& track); 0207 0208 //________________________________________________________________ 0209 // Process information 0210 struct G4ITBrownianState : public G4ITTransportationState 0211 { 0212 public: 0213 G4ITBrownianState(); 0214 ~G4ITBrownianState() override 0215 { 0216 ; 0217 } 0218 G4String GetType() override 0219 { 0220 return "G4ITBrownianState"; 0221 } 0222 0223 G4bool fPathLengthWasCorrected; 0224 G4bool fTimeStepReachedLimit; 0225 G4bool fComputeLastPosition; 0226 G4double fRandomNumber; 0227 }; 0228 0229 G4bool fUseMaximumTimeBeforeReachingBoundary; 0230 G4Material* fNistWater; 0231 0232 G4bool fUseSchedulerMinTimeSteps; 0233 G4double fInternalMinTimeStep; 0234 G4bool fSpeedMeUp; 0235 0236 // Water density table 0237 const std::vector<G4double>* fpWaterDensity; 0238 0239 G4BrownianAction* fpBrownianAction; 0240 G4VUserBrownianAction *fpUserBrownianAction; 0241 0242 }; 0243 0244 0245 inline void G4DNABrownianTransportation::SetBrownianAction(G4BrownianAction* brownianAction) 0246 { 0247 fpBrownianAction = brownianAction; 0248 } 0249 0250 inline void G4DNABrownianTransportation::SetUserBrownianAction(G4VUserBrownianAction* brownianAction) 0251 { 0252 fpUserBrownianAction = brownianAction; 0253 } 0254 0255 0256 #endif // G4ITBROWNIANTRANSPORTATION_H
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |