Back to home page

EIC code displayed by LXR

 
 

    


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

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 /// file: G4DNAPolyNucleotideReactionProcess.hh
0027 /// brief: This file handls reaction process with DNA geometry.
0028 #ifndef G4DNAPolyNucleotideReactionProcess_hh
0029 #define G4DNAPolyNucleotideReactionProcess_hh
0030 #include <CLHEP/Units/SystemOfUnits.h>
0031 #include "G4VITDiscreteProcess.hh"
0032 #include <variant>
0033 class G4DNAMolecularReactionTable;
0034 class G4DNAComponentNode;
0035 class G4VDNAHitModel;
0036 
0037 class G4DNAPolyNucleotideReactionProcess : public G4VITDiscreteProcess
0038 {
0039  public:
0040   explicit G4DNAPolyNucleotideReactionProcess(
0041     const G4String& aName = "DNAStaticMoleculeReactionProcess",
0042     G4int verbosityLevel  = 0);
0043   ~G4DNAPolyNucleotideReactionProcess() override;
0044   
0045   G4DNAPolyNucleotideReactionProcess(const G4DNAPolyNucleotideReactionProcess&) = delete;
0046   G4DNAPolyNucleotideReactionProcess& operator =(
0047     const G4DNAPolyNucleotideReactionProcess&) = delete;
0048  
0049   inline void SetDNADamageReactionModel(G4VDNAHitModel* pModel);
0050 
0051   G4bool IsApplicable(const G4ParticleDefinition&) override { return true; }
0052 
0053   G4double CalculateTimeStep(const G4Track& trackA,
0054                              const G4double& userTimeStep = 0);
0055 
0056   void StartTracking(G4Track* aTrack) override;
0057 
0058   G4double PostStepGetPhysicalInteractionLength(
0059     const G4Track&,  // track
0060     G4double,        // previousStepSize
0061     G4ForceCondition* pForceCond) override;
0062   G4VParticleChange* PostStepDoIt(const G4Track& track, const G4Step&) override;
0063 
0064   G4double GetMeanFreePath(const G4Track&, G4double, G4ForceCondition*) override
0065   {
0066     return DBL_MAX;
0067   }
0068 
0069   inline void SetVerbose(G4int verbose);
0070 
0071  protected:
0072   G4VParticleChange fParticleChange;
0073 
0074  private:
0075   struct G4PolyNucleotideReactionState
0076     : public G4ProcessStateBase<G4PolyNucleotideReactionState>
0077   {
0078     G4PolyNucleotideReactionState();
0079     ~G4PolyNucleotideReactionState() override = default;
0080     G4String GetType() override { return "PolyNucleotideReactionState"; }
0081 
0082     using DNANode =
0083       std::variant<const G4DNAComponentNode*, /*for dnadamage chain*/
0084                    const G4VPhysicalVolume* /*for molecularDNA chain*/>;
0085     DNANode fNodeReactant;
0086     G4double fSampledMinTimeStep;
0087     G4double fPreviousTimeAtPreStepPoint;
0088   };
0089   G4bool fHasAlreadyReachedNullTime{false};
0090   G4int fVerbose;
0091   G4double fRCutOff;
0092   G4VDNAHitModel* fpDamageModel{nullptr};
0093 };
0094 inline void G4DNAPolyNucleotideReactionProcess::SetVerbose(G4int verbose)
0095 {
0096   fVerbose = verbose;
0097 }
0098 
0099 inline void G4DNAPolyNucleotideReactionProcess::SetDNADamageReactionModel(
0100   G4VDNAHitModel* pModel)
0101 {
0102   fpDamageModel = pModel;
0103 }
0104 #endif