Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAMolecularDissociation Class Reference

#include <G4DNAMolecularDissociation.hh>

+ Inheritance diagram for G4DNAMolecularDissociation:

Public Types

using Species = const G4MoleculeDefinition
 
using Displacer = G4VMolecularDissociationDisplacer
 

Public Member Functions

 G4DNAMolecularDissociation (const G4String &processName="DNAMolecularDecay", G4ProcessType type=fDecay)
 
 G4DNAMolecularDissociation ()=delete
 
 G4DNAMolecularDissociation (const G4DNAMolecularDissociation &right)=delete
 
G4DNAMolecularDissociationoperator= (const G4DNAMolecularDissociation &right)=delete
 
 ~G4DNAMolecularDissociation () override
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition) override
 
G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &step) override
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &step) override
 
void SetVerbose (G4int)
 
void SetUserBrownianAction (G4VUserBrownianAction *pBrownianAction)
 
void SetDisplacer (Species *, Displacer *)
 
DisplacerGetDisplacer (Species *)
 
- Public Member Functions inherited from G4VITRestDiscreteProcess
 G4VITRestDiscreteProcess ()=delete
 
 G4VITRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VITRestDiscreteProcess (const G4VITRestDiscreteProcess &)=delete
 
G4VITRestDiscreteProcessoperator= (const G4VITRestDiscreteProcess &right)=delete
 
 ~G4VITRestDiscreteProcess () override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *) override
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
- Public Member Functions inherited from G4VITProcess
 G4VITProcess (const G4String &name, G4ProcessType type=fNotDefined)
 
virtual ~G4VITProcess ()
 
 G4VITProcess (const G4VITProcess &other)
 
G4VITProcessoperator= (const G4VITProcess &other)
 
G4bool operator== (const G4VITProcess &right) const
 
G4bool operator!= (const G4VITProcess &right) const
 
size_t GetProcessID () const
 
G4shared_ptr< G4ProcessState_LockGetProcessState ()
 
void SetProcessState (G4shared_ptr< G4ProcessState_Lock > aProcInfo)
 
void ResetProcessState ()
 
virtual void StartTracking (G4Track *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double GetInteractionTimeLeft ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4bool ProposesTimeStep () const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

virtual G4VParticleChangeDecayIt (const G4Track &, const G4Step &)
 
G4double GetMeanLifeTime (const G4Track &, G4ForceCondition *) override
 
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *) override
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VITProcess
void RetrieveProcessInfo ()
 
void CreateInfo ()
 
template<typename T >
T * GetState ()
 
virtual void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
virtual void ClearInteractionTimeLeft ()
 
virtual void ClearNumberOfInteractionLengthLeft ()
 
void SetInstantiateProcessState (G4bool flag)
 
G4bool InstantiateProcessState ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VITProcess
static const size_t & GetMaxProcessIndex ()
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VITProcess
G4shared_ptr< G4ProcessStatefpState
 
G4bool fProposesTimeStep
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

G4DNAMolecularDissociation should be called only for molecules. It will dissociate the molecules using the decay associated to this molecule and if a displacement scheme has been registered, it will place the products to the expected position.

Definition at line 63 of file G4DNAMolecularDissociation.hh.

Member Typedef Documentation

◆ Displacer

◆ Species

Constructor & Destructor Documentation

◆ G4DNAMolecularDissociation() [1/3]

G4DNAMolecularDissociation::G4DNAMolecularDissociation ( const G4String processName = "DNAMolecularDecay",
G4ProcessType  type = fDecay 
)

Definition at line 50 of file G4DNAMolecularDissociation.cc.

53 : G4VITRestDiscreteProcess(processName, type)
54{
55 // set Process Sub Type
57 enableAlongStepDoIt = false;
58 enablePostStepDoIt = true;
59 enableAtRestDoIt = true;
60
61 fVerbose = 0;
62
63#ifdef G4VERBOSE
64 if (verboseLevel > 1)
65 {
66 G4cout << "G4MolecularDissociationProcess constructor " << " Name:"
67 << processName << G4endl;
68 }
69#endif
70
72
73 fDecayAtFixedTime = true;
74 fProposesTimeStep = true;
75}
@ fLowEnergyMolecularDecay
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4bool fProposesTimeStep
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:363
G4int verboseLevel
Definition: G4VProcess.hh:360
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:364
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:365
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325

◆ G4DNAMolecularDissociation() [2/3]

G4DNAMolecularDissociation::G4DNAMolecularDissociation ( )
delete

◆ G4DNAMolecularDissociation() [3/3]

G4DNAMolecularDissociation::G4DNAMolecularDissociation ( const G4DNAMolecularDissociation right)
delete

◆ ~G4DNAMolecularDissociation()

G4DNAMolecularDissociation::~G4DNAMolecularDissociation ( )
override

Definition at line 79 of file G4DNAMolecularDissociation.cc.

80{
81 delete fpBrownianAction;
82}

Member Function Documentation

◆ AtRestDoIt()

G4VParticleChange * G4DNAMolecularDissociation::AtRestDoIt ( const G4Track track,
const G4Step step 
)
overridevirtual

Reimplemented from G4VITRestDiscreteProcess.

Definition at line 336 of file G4DNAMolecularDissociation.cc.

338{
341 return DecayIt(track, step);
342}
virtual G4VParticleChange * DecayIt(const G4Track &, const G4Step &)
virtual void ClearInteractionTimeLeft()
virtual void ClearNumberOfInteractionLengthLeft()

Referenced by PostStepDoIt().

◆ AtRestGetPhysicalInteractionLength()

G4double G4DNAMolecularDissociation::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
overridevirtual

Reimplemented from G4VITRestDiscreteProcess.

Definition at line 346 of file G4DNAMolecularDissociation.cc.

348{
349 if (fDecayAtFixedTime)
350 {
351 return GetMeanLifeTime(track, condition);
352 }
353
355}
G4double condition(const G4ErrorSymMatrix &m)
G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *) override
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *) override

◆ DecayIt()

G4VParticleChange * G4DNAMolecularDissociation::DecayIt ( const G4Track track,
const G4Step  
)
protectedvirtual

Definition at line 119 of file G4DNAMolecularDissociation.cc.

121{
123 auto pMotherMolecule = GetMolecule(track);
124 auto pMotherMoleculeDefinition = pMotherMolecule->GetDefinition();
125
126 if (pMotherMoleculeDefinition->GetDecayTable())
127 {
128 const auto pDissociationChannels = pMotherMolecule->GetDissociationChannels();
129
130 if (pDissociationChannels == nullptr)
131 {
132 G4ExceptionDescription exceptionDescription;
133 pMotherMolecule->PrintState();
134 exceptionDescription << "No decay channel was found for the molecule : "
135 << pMotherMolecule->GetName() << G4endl;
136 G4Exception("G4DNAMolecularDissociation::DecayIt",
137 "G4DNAMolecularDissociation::NoDecayChannel",
139 exceptionDescription);
140 return &aParticleChange;
141 }
142
143 auto decayVectorSize = pDissociationChannels->size();
144 G4double RdmValue = G4UniformRand();
145
146 const G4MolecularDissociationChannel* pDecayChannel = nullptr;
147 size_t i = 0;
148 do
149 {
150 pDecayChannel = (*pDissociationChannels)[i];
151 if (RdmValue < pDecayChannel->GetProbability())
152 {
153 break;
154 }
155 RdmValue -= pDecayChannel->GetProbability();
156 i++;
157 } while (i < decayVectorSize);
158
159 G4double decayEnergy = pDecayChannel->GetEnergy();
160 auto nbProducts = pDecayChannel->GetNbProducts();
161
162 if (decayEnergy > 0.)
163 {
165 }
166
167 if (nbProducts)
168 {
169 std::vector<G4ThreeVector> productsDisplacement(nbProducts);
170 G4ThreeVector motherMoleculeDisplacement;
171
172 auto it = fDisplacementMap.find(pMotherMoleculeDefinition);
173
174 if (it != fDisplacementMap.end())
175 {
176 auto pDisplacer = it->second.get();
177 productsDisplacement = pDisplacer->GetProductsDisplacement(pDecayChannel);
178 motherMoleculeDisplacement =
179 pDisplacer->GetMotherMoleculeDisplacement(pDecayChannel);
180 }
181 else
182 {
184 errMsg << "No G4MolecularDecayProcess::theDecayDisplacementMap["
185 << pMotherMolecule->GetName() + "]";
186 G4Exception("G4MolecularDecayProcess::DecayIt",
187 "DNAMolecularDecay001",
189 errMsg);
190 }
191
193
194#ifdef G4VERBOSE
195 if (fVerbose)
196 {
197 G4cout << "Decay Process : " << pMotherMolecule->GetName()
198 << " (trackID :" << track.GetTrackID() << ") "
199 << pDecayChannel->GetName() << G4endl;
200 }
201#endif
202
204
205 for (G4int j = 0; j < nbProducts; j++)
206 {
207 auto pProduct = new G4Molecule(pDecayChannel->GetProduct(j));
208
209 G4ThreeVector displacement = motherMoleculeDisplacement + productsDisplacement[j];
210 double mag_displacement = displacement.mag();
211 G4ThreeVector displacement_direction = displacement / (mag_displacement + 1e-30);
212
213 double prNewSafety = DBL_MAX;
214
215 //double step =
216 pNavigator->CheckNextStep(track.GetPosition(),
217 displacement_direction,
218 mag_displacement,
219 prNewSafety);
220
221 //if(prNewSafety < mag_displacement || step < mag_displacement)
222 mag_displacement = std::min(prNewSafety * 0.8, mag_displacement);
223
224 G4ThreeVector product_pos = track.GetPosition()
225 + displacement_direction * mag_displacement;
226
227 //Hoang: force changing position track::
228 if(fpBrownianAction != nullptr)
229 {
230 fpBrownianAction->Transport(product_pos);
231 }
232 //Hoang: force changing position track
233
234 const G4AffineTransform& transform = pNavigator->GetGlobalToLocalTransform();
235
236 G4ThreeVector localPoint = transform.TransformPoint(product_pos); //track.GetPosition());
237 // warning if the decayed product is outside of the volume and
238 // the mother volume has no water material (the decayed product
239 // is outside of the world volume will be killed in the next step)
240 if (track.GetTouchable()->GetSolid()->Inside(localPoint) !=
241 EInside::kInside)
242 {
243 auto WaterMaterial = G4Material::GetMaterial("G4_WATER");
244 auto Motherlogic = track.GetTouchable()->GetVolume()->
245 GetMotherLogical();
246 if (Motherlogic != nullptr
247 && Motherlogic->GetMaterial() != WaterMaterial)
248 {
250 ED << "The decayed product is outside of the volume : "
251 << track.GetTouchable()->GetVolume()->GetName()
252 << " with material : "<< Motherlogic->GetMaterial()
253 ->GetName()<< G4endl;
254 G4Exception("G4DNAMolecularDissociation::DecayIt()",
255 "OUTSIDE_OF_MOTHER_VOLUME",
256 JustWarning, ED);
257 }
258 }
259
260 auto pSecondary = pProduct->BuildTrack(track.GetGlobalTime(), product_pos);
261
262 pSecondary->SetTrackStatus(fAlive);
263#ifdef G4VERBOSE
264 if (fVerbose)
265 {
266 G4cout << "Product : " << pProduct->GetName() << G4endl;
267 }
268#endif
269 // add the secondary track in the List
270 aParticleChange.G4VParticleChange::AddSecondary(pSecondary);
271 }
272#ifdef G4VERBOSE
273 if (fVerbose)
274 {
275 G4cout << "-------------" << G4endl;
276 }
277#endif
278 }
279 //DEBUG
280 else if (fVerbose && decayEnergy)
281 {
282 G4cout << "No products for this channel" << G4endl;
283 G4cout << "-------------" << G4endl;
284 }
285 /*
286 else if(!decayEnergy && !nbProducts)
287 {
288 G4ExceptionDescription errMsg;
289 errMsg << "There is no products and no energy specified in the molecular "
290 "dissociation channel";
291 G4Exception("G4MolecularDissociationProcess::DecayIt",
292 "DNAMolecularDissociation002",
293 FatalErrorInArgument,
294 errMsg);
295 }
296 */
297 }
298
300
301 return &aParticleChange;
302}
@ JustWarning
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:74
@ fAlive
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
double mag() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
static G4ITTransportationManager * GetTransportationManager()
G4ITNavigator * GetNavigatorForTracking() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:691
void Initialize(const G4Track &) override
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
const G4VTouchable * GetTouchable() const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual G4VSolid * GetSolid(G4int depth=0) const
Definition: G4VTouchable.cc:42
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34
virtual void Transport(G4ThreeVector &, G4Track *pTrack=nullptr)=0
#define DBL_MAX
Definition: templates.hh:62

Referenced by AtRestDoIt().

◆ GetDisplacer()

G4DNAMolecularDissociation::Displacer * G4DNAMolecularDissociation::GetDisplacer ( Species pSpecies)

Definition at line 313 of file G4DNAMolecularDissociation.cc.

314{
315 return fDisplacementMap[pSpecies].get();
316}

◆ GetMeanFreePath()

G4double G4DNAMolecularDissociation::GetMeanFreePath ( const G4Track ,
G4double  ,
G4ForceCondition  
)
overrideprotectedvirtual

Implements G4VITRestDiscreteProcess.

Definition at line 367 of file G4DNAMolecularDissociation.cc.

370{
371 return 0;
372}

◆ GetMeanLifeTime()

G4double G4DNAMolecularDissociation::GetMeanLifeTime ( const G4Track track,
G4ForceCondition  
)
overrideprotectedvirtual

Implements G4VITRestDiscreteProcess.

Definition at line 110 of file G4DNAMolecularDissociation.cc.

112{
113 G4double output = GetMolecule(track)->GetDecayTime() - track.GetProperTime();
114 return output > 0. ? output : 0.;
115}
G4double GetDecayTime() const
Definition: G4Molecule.cc:472
G4double GetProperTime() const

Referenced by AtRestGetPhysicalInteractionLength().

◆ IsApplicable()

G4bool G4DNAMolecularDissociation::IsApplicable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 86 of file G4DNAMolecularDissociation.cc.

88{
89 if (aParticleType.GetParticleType() == "Molecule")
90 {
91#ifdef G4VERBOSE
92
93 if (fVerbose > 1)
94 {
95 G4cout << "G4MolecularDissociation::IsApplicable(";
96 G4cout << aParticleType.GetParticleName() << ",";
97 G4cout << aParticleType.GetParticleType() << ")" << G4endl;
98 }
99#endif
100 return (true);
101 }
102 else
103 {
104 return false;
105 }
106}
const G4String & GetParticleType() const
const G4String & GetParticleName() const

◆ operator=()

G4DNAMolecularDissociation & G4DNAMolecularDissociation::operator= ( const G4DNAMolecularDissociation right)
delete

◆ PostStepDoIt()

G4VParticleChange * G4DNAMolecularDissociation::PostStepDoIt ( const G4Track track,
const G4Step step 
)
overridevirtual

Reimplemented from G4VITRestDiscreteProcess.

Definition at line 359 of file G4DNAMolecularDissociation.cc.

361{
362 return AtRestDoIt(track, step);
363}
G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &step) override

◆ PostStepGetPhysicalInteractionLength()

G4double G4DNAMolecularDissociation::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
overridevirtual

Reimplemented from G4VITRestDiscreteProcess.

Definition at line 320 of file G4DNAMolecularDissociation.cc.

323{
324 return 0; //-1*(-log(1-G4UniformRand())*100*1e-15*s);
325}

◆ SetDisplacer()

void G4DNAMolecularDissociation::SetDisplacer ( Species pSpecies,
Displacer pDisplacer 
)

Definition at line 306 of file G4DNAMolecularDissociation.cc.

307{
308 fDisplacementMap.emplace(pSpecies, std::unique_ptr<Displacer>(pDisplacer));
309}

Referenced by G4EmDNAChemistry::ConstructProcess(), G4EmDNAChemistry_option1::ConstructProcess(), and G4EmDNAChemistry_option2::ConstructProcess().

◆ SetUserBrownianAction()

void G4DNAMolecularDissociation::SetUserBrownianAction ( G4VUserBrownianAction pBrownianAction)
inline

Definition at line 91 of file G4DNAMolecularDissociation.hh.

92 {
93 fpBrownianAction = pBrownianAction;
94 }

◆ SetVerbose()

void G4DNAMolecularDissociation::SetVerbose ( G4int  verbose)

Definition at line 329 of file G4DNAMolecularDissociation.cc.

330{
331 fVerbose = verbose;
332}

The documentation for this class was generated from the following files: