Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAMoleculeEncounterStepper.hh
Go to the documentation of this file.
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// $Id: G4DNAMoleculeEncounterStepper.hh 64057 2012-10-30 15:04:49Z gcosmo $
27//
28// Author: Mathieu Karamitros ([email protected])
29//
30// WARNING : This class is released as a prototype.
31// It might strongly evolve or even disapear in the next releases.
32//
33// History:
34// -----------
35// 10 Oct 2011 M.Karamitros created
36//
37// -------------------------------------------------------------------
38
39#ifndef G4MOLECULEENCOUNTERSTEPPER_H
40#define G4MOLECULEENCOUNTERSTEPPER_H
41
42#include "G4VITTimeStepper.hh"
43#include "G4ITManager.hh"
44
47
48class G4Molecule;
49
50/**
51 * Given a molecule G4DNAMoleculeEncounterStepper will calculate for its possible reactants
52 * what will be the minimum encounter time and the associated molecules.*
53 *
54 * This model includes dynamical time steps as explained in
55 * "Computer-Aided Stochastic Modeling of the Radiolysis of Liquid Water",
56 * V. Michalik, M. Begusová, E. A. Bigildeev,
57 * Radiation Research, Vol. 149, No. 3 (Mar., 1998), pp. 224-236
58 *
59 */
60
62{
63public:
68
69 virtual void Prepare();
70// virtual void PrepareForAllProcessors();
71 virtual G4double CalculateStep(const G4Track&, const G4double&);
72
75
76 inline void SetVerbose(int);
77 // Final time returned when reaction is avalaible in the reaction table = 1
78 // All details = 2
79
80private:
81
82 void RetrieveResults(const G4Track&, const G4Molecule*, const G4Molecule*, const G4double /*reactionRange*/,
83 G4KDTreeResultHandle&, G4bool iterate = true);
84
85 G4bool fHasAlreadyReachedNullTime;
86
88 const G4DNAMolecularReactionTable*& fMolecularReactionTable ;
89 G4VDNAReactionModel* fReactionModel;
90 G4int fVerbose ;
91};
92
94{
95 fReactionModel = reactionModel ;
96}
97
99{
100 return fReactionModel;
101}
102
104{
105 fVerbose = flag;
106}
107
108#endif // G4MOLECULEENCOUNTERSTEPPER_H
#define G4IT_ADD_CLONE(parent_class, kid_class)
Definition: AddClone_def.hh:45
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
virtual G4double CalculateStep(const G4Track &, const G4double &)
void SetReactionModel(G4VDNAReactionModel *)