Geant4  10.02.p03
G4EmDNAChemistry.cc
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 #include "G4EmDNAChemistry.hh"
27 
28 #include "G4PhysicalConstants.hh"
29 #include "G4SystemOfUnits.hh"
30 
32 #include "G4DNAChemistryManager.hh"
33 #include "G4DNAWaterExcitationStructure.hh"
34 #include "G4ProcessManager.hh"
35 
37 
38 // *** Processes and models for Geant4-DNA
39 
41 
42 #include "G4DNAAttachment.hh"
43 #include "G4DNAVibExcitation.hh"
44 
45 #include "G4DNAElastic.hh"
49 
54 #include "G4VDNAReactionModel.hh"
56 
58 
59 // particles
60 
61 #include "G4Electron.hh"
62 #include "G4Proton.hh"
63 #include "G4GenericIon.hh"
64 
65 #include "G4MoleculeTable.hh"
66 #include "G4H2O.hh"
67 #include "G4H2.hh"
68 #include "G4Hydrogen.hh"
69 #include "G4OH.hh"
70 #include "G4H3O.hh"
71 #include "G4Electron_aq.hh"
72 #include "G4H2O2.hh"
73 
74 #include "G4PhysicsListHelper.hh"
75 #include "G4BuilderType.hh"
76 
77 /****/
79 #include "G4ProcessVector.hh"
80 #include "G4ProcessTable.hh"
83 /****/
84 
85 // factory
87 
89 
90 #include "G4Threading.hh"
91 
94 {
96 }
97 
98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
99 
101 {
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105 
107 {
108  //-----------------------------------
109  G4Electron::Definition(); // safety
110 
111  //-----------------------------------
112  // Create the definition
120 
121 // G4MoleculeTable::Instance()->CreateMoleculeModel("H3Op", G4H3O::Definition());
122 // G4Molecule* OHm = G4MoleculeTable::Instance()->
123 // CreateMoleculeModel("OHm", // just a tag to store and retrieve from
124 // // G4MoleculeTable
125 // G4OH::Definition(), -1, 5.0e-9 * (m2 / s));
126 // OHm->SetMass(17.0079 * g / Avogadro * c_squared);
127 // G4MoleculeTable::Instance()->CreateMoleculeModel("OH", G4OH::Definition());
128 // G4MoleculeTable::Instance()->CreateMoleculeModel("e_aq",
129 // G4Electron_aq::Definition());
130 // G4MoleculeTable::Instance()->CreateMoleculeModel("H",
131 // G4Hydrogen::Definition());
132 // G4MoleculeTable::Instance()->CreateMoleculeModel("H2", G4H2::Definition());
133 // G4MoleculeTable::Instance()->CreateMoleculeModel("H2O2", G4H2O2::Definition());
134 
135  //____________________________________________________________________________
136 
139  CreateConfiguration("OHm", // just a tag to store and retrieve from
140  // G4MoleculeTable
142  -1, // charge
143  5.0e-9 * (m2 / s));
144  OHm->SetMass(17.0079 * g / Avogadro * c_squared);
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155 
157 {
158 // //-----------------------------------
159 // //Create the dynamic objects
160 // G4Molecule* OH = G4MoleculeTable::Instance()->GetMoleculeModel("OH");
161 // G4Molecule* OHm = G4MoleculeTable::Instance()->GetMoleculeModel("OHm");
162 // G4Molecule* e_aq = G4MoleculeTable::Instance()->GetMoleculeModel("e_aq");
163 // G4Molecule* H2 = G4MoleculeTable::Instance()->GetMoleculeModel("H2");
164 // G4Molecule* H3O = G4MoleculeTable::Instance()->GetMoleculeModel("H3Op");
165 // G4Molecule* H = G4MoleculeTable::Instance()->GetMoleculeModel("H");
166  //-----------------------------------
167  //Get the molecular configuration
180 
181  //-------------------------------------
182  //Define the decay channels
186 
188  *(water->GetGroundStateElectronOccupancy()));
189 
191  // EXCITATIONS //
193  G4DNAWaterExcitationStructure waterExcitation;
194  //--------------------------------------------------------
195  //---------------Excitation on the fifth layer------------
196 
197  decCh1 = new G4MolecularDissociationChannel("A^1B_1_Relaxation");
198  decCh2 = new G4MolecularDissociationChannel("A^1B_1_DissociativeDecay");
199  //Decay 1 : OH + H
200  decCh1->SetEnergy(waterExcitation.ExcitationEnergy(0));
201  decCh1->SetProbability(0.35);
203 
204  decCh2->AddProduct(OH);
205  decCh2->AddProduct(H);
206  decCh2->SetProbability(0.65);
207  decCh2->SetDisplacementType(
209 
210 // water->AddExcitedState("A^1B_1");
211  occ->RemoveElectron(4, 1); // this is the transition form ground state to
212  occ->AddElectron(5, 1); // the first unoccupied orbital: A^1B_1
213 
214  water->NewConfigurationWithElectronOccupancy("A^1B_1", *occ);
215  water->AddDecayChannel("A^1B_1", decCh1);
216  water->AddDecayChannel("A^1B_1", decCh2);
217 
218  //--------------------------------------------------------
219  //---------------Excitation on the fourth layer-----------
220  decCh1 = new G4MolecularDissociationChannel("B^1A_1_Relaxation_Channel");
221  decCh2 = new G4MolecularDissociationChannel("B^1A_1_DissociativeDecay");
223  "B^1A_1_AutoIonisation_Channel");
224 
225  //Decay 1 : energy
226  decCh1->SetEnergy(waterExcitation.ExcitationEnergy(1));
227  decCh1->SetProbability(0.3);
228 
229  //Decay 2 : 2OH + H_2
230  decCh2->AddProduct(H2);
231  decCh2->AddProduct(OH);
232  decCh2->AddProduct(OH);
233  decCh2->SetProbability(0.15);
234  decCh2->SetDisplacementType(
236 
237  //Decay 3 : OH + H_3Op + e_aq
238  decCh3->AddProduct(OH);
239  decCh3->AddProduct(H3O);
240  decCh3->AddProduct(e_aq);
241  decCh3->SetProbability(0.55);
243 
244  *occ = *(water->GetGroundStateElectronOccupancy());
245  occ->RemoveElectron(3); // this is the transition form ground state to
246  occ->AddElectron(5, 1); // the first unoccupied orbital: B^1A_1
247 
248  water->NewConfigurationWithElectronOccupancy("B^1A_1", *occ);
249  water->AddDecayChannel("B^1A_1", decCh1);
250  water->AddDecayChannel("B^1A_1", decCh2);
251  water->AddDecayChannel("B^1A_1", decCh3);
252 
253  //-------------------------------------------------------
254  //-------------------Excitation of 3rd layer-----------------
255  decCh1 = new G4MolecularDissociationChannel(
256  "Excitation3rdLayer_AutoIonisation_Channel");
257  decCh2 = new G4MolecularDissociationChannel(
258  "Excitation3rdLayer_Relaxation_Channel");
259 
260  //Decay channel 1 : : OH + H_3Op + e_aq
261  decCh1->AddProduct(OH);
262  decCh1->AddProduct(H3O);
263  decCh1->AddProduct(e_aq);
264 
265  decCh1->SetProbability(0.5);
267 
268  //Decay channel 2 : energy
269  decCh2->SetEnergy(waterExcitation.ExcitationEnergy(2));
270  decCh2->SetProbability(0.5);
271 
272  //Electronic configuration of this decay
273  *occ = *(water->GetGroundStateElectronOccupancy());
274  occ->RemoveElectron(2, 1);
275  occ->AddElectron(5, 1);
276 
277  //Configure the water molecule
278  water->NewConfigurationWithElectronOccupancy("Excitation3rdLayer", *occ);
279  water->AddDecayChannel("Excitation3rdLayer", decCh1);
280  water->AddDecayChannel("Excitation3rdLayer", decCh2);
281 
282  //-------------------------------------------------------
283  //-------------------Excitation of 2nd layer-----------------
284  decCh1 = new G4MolecularDissociationChannel(
285  "Excitation2ndLayer_AutoIonisation_Channel");
286  decCh2 = new G4MolecularDissociationChannel(
287  "Excitation2ndLayer_Relaxation_Channel");
288 
289  //Decay Channel 1 : : OH + H_3Op + e_aq
290  decCh1->AddProduct(OH);
291  decCh1->AddProduct(H3O);
292  decCh1->AddProduct(e_aq);
293 
294  decCh1->SetProbability(0.5);
296 
297  //Decay channel 2 : energy
298  decCh2->SetEnergy(waterExcitation.ExcitationEnergy(3));
299  decCh2->SetProbability(0.5);
300 
301  *occ = *(water->GetGroundStateElectronOccupancy());
302  occ->RemoveElectron(1, 1);
303  occ->AddElectron(5, 1);
304 
305  water->NewConfigurationWithElectronOccupancy("Excitation2ndLayer", *occ);
306  water->AddDecayChannel("Excitation2ndLayer", decCh1);
307  water->AddDecayChannel("Excitation2ndLayer", decCh2);
308 
309  //-------------------------------------------------------
310  //-------------------Excitation of 1st layer-----------------
311  decCh1 = new G4MolecularDissociationChannel(
312  "Excitation1stLayer_AutoIonisation_Channel");
313  decCh2 = new G4MolecularDissociationChannel(
314  "Excitation1stLayer_Relaxation_Channel");
315 
316  *occ = *(water->GetGroundStateElectronOccupancy());
317  occ->RemoveElectron(0, 1);
318  occ->AddElectron(5, 1);
319 
320  //Decay Channel 1 : : OH + H_3Op + e_aq
321  decCh1->AddProduct(OH);
322  decCh1->AddProduct(H3O);
323  decCh1->AddProduct(e_aq);
324  decCh1->SetProbability(0.5);
326 
327  //Decay channel 2 : energy
328  decCh2->SetEnergy(waterExcitation.ExcitationEnergy(4));
329  decCh2->SetProbability(0.5);
330 
331  water->NewConfigurationWithElectronOccupancy("Excitation1stLayer", *occ);
332  water->AddDecayChannel("Excitation1stLayer", decCh1);
333  water->AddDecayChannel("Excitation1stLayer", decCh2);
334 
336  // IONISATION //
338  //--------------------------------------------------------
339  //------------------- Ionisation -------------------------
340 
341  decCh1 = new G4MolecularDissociationChannel("Ionisation_Channel");
342 
343  //Decay Channel 1 : : OH + H_3Op
344  decCh1->AddProduct(H3O);
345  decCh1->AddProduct(OH);
346  decCh1->SetProbability(1);
347  decCh1->SetDisplacementType(
349 
350  *occ = *(water->GetGroundStateElectronOccupancy());
351  occ->RemoveElectron(4, 1);
352  // this is a ionized h2O with a hole in its last orbital
353  water->NewConfigurationWithElectronOccupancy("Ionisation5", *occ);
354  water->AddDecayChannel("Ionisation5",
355  decCh1);
356 
357  *occ = *(water->GetGroundStateElectronOccupancy());
358  occ->RemoveElectron(3, 1);
359  water->NewConfigurationWithElectronOccupancy("Ionisation4", *occ);
360  water->AddDecayChannel("Ionisation4",
361  new G4MolecularDissociationChannel(*decCh1));
362 
363  *occ = *(water->GetGroundStateElectronOccupancy());
364  occ->RemoveElectron(2, 1);
365  water->NewConfigurationWithElectronOccupancy("Ionisation3", *occ);
366  water->AddDecayChannel("Ionisation3",
367  new G4MolecularDissociationChannel(*decCh1));
368 
369  *occ = *(water->GetGroundStateElectronOccupancy());
370  occ->RemoveElectron(1, 1);
371  water->NewConfigurationWithElectronOccupancy("Ionisation2", *occ);
372  water->AddDecayChannel("Ionisation2",
373  new G4MolecularDissociationChannel(*decCh1));
374 
375  *occ = *(water->GetGroundStateElectronOccupancy());
376  occ->RemoveElectron(0, 1);
377  water->NewConfigurationWithElectronOccupancy("Ionisation1", *occ);
378  water->AddDecayChannel("Ionisation1",
379  new G4MolecularDissociationChannel(*decCh1));
380 
382  // Dissociative Attachment //
384  decCh1 = new G4MolecularDissociationChannel("DissociativeAttachment");
385 
386  //Decay 1 : 2OH + H_2
387  decCh1->AddProduct(H2);
388  decCh1->AddProduct(OHm);
389  decCh1->AddProduct(OH);
390  decCh1->SetProbability(1);
392  DissociativeAttachment);
393 
394  *occ = *(water->GetGroundStateElectronOccupancy());
395  occ->AddElectron(5, 1); // H_2O^-
396  water->NewConfigurationWithElectronOccupancy("DissociativeAttachment", *occ);
397  water->AddDecayChannel("DissociativeAttachment", decCh1);
398 
399  delete occ;
400 }
401 
402 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
403 
405  theReactionTable)
406 {
407  // //-----------------------------------
408 // G4Molecule* OH = G4MoleculeTable::Instance()->GetMoleculeModel("OH");
409 // G4Molecule* e_aq = G4MoleculeTable::Instance()->GetMoleculeModel("e_aq");
410 // G4Molecule* H2 = G4MoleculeTable::Instance()->GetMoleculeModel("H2");
411 // G4Molecule* H3Op = G4MoleculeTable::Instance()->GetMoleculeModel("H3Op");
412 // G4Molecule* OHm = G4MoleculeTable::Instance()->GetMoleculeModel("OHm");
413 // G4Molecule* H2O2 = G4MoleculeTable::Instance()->GetMoleculeModel("H2O2");
414 // G4Molecule* H = G4MoleculeTable::Instance()->GetMoleculeModel("H");
415  //-----------------------------------
416  //Get the molecular configuration
431 
432  //------------------------------------------------------------------
433  // e_aq + e_aq + 2H2O -> H2 + 2OH-
435  0.5e10 * (1e-3 * m3 / (mole * s)), e_aq, e_aq);
436  reactionData->AddProduct(OHm);
437  reactionData->AddProduct(OHm);
438  reactionData->AddProduct(H2);
439  theReactionTable->SetReaction(reactionData);
440  //------------------------------------------------------------------
441  // e_aq + *OH -> OH-
442  reactionData = new G4DNAMolecularReactionData(
443  2.95e10 * (1e-3 * m3 / (mole * s)), e_aq, OH);
444  reactionData->AddProduct(OHm);
445  theReactionTable->SetReaction(reactionData);
446  //------------------------------------------------------------------
447  // e_aq + H* + H2O -> H2 + OH-
448  reactionData = new G4DNAMolecularReactionData(
449  2.65e10 * (1e-3 * m3 / (mole * s)), e_aq, H);
450  reactionData->AddProduct(OHm);
451  reactionData->AddProduct(H2);
452  theReactionTable->SetReaction(reactionData);
453  //------------------------------------------------------------------
454  // e_aq + H3O+ -> H* + H2O
455  reactionData = new G4DNAMolecularReactionData(
456  2.11e10 * (1e-3 * m3 / (mole * s)), e_aq, H3Op);
457  reactionData->AddProduct(H);
458  theReactionTable->SetReaction(reactionData);
459  //------------------------------------------------------------------
460  // e_aq + H2O2 -> OH- + *OH
461  reactionData = new G4DNAMolecularReactionData(
462  1.41e10 * (1e-3 * m3 / (mole * s)), e_aq, H2O2);
463  reactionData->AddProduct(OHm);
464  reactionData->AddProduct(OH);
465  theReactionTable->SetReaction(reactionData);
466  //------------------------------------------------------------------
467  // *OH + *OH -> H2O2
468  reactionData = new G4DNAMolecularReactionData(
469  0.44e10 * (1e-3 * m3 / (mole * s)), OH, OH);
470  reactionData->AddProduct(H2O2);
471  theReactionTable->SetReaction(reactionData);
472  //------------------------------------------------------------------
473  // *OH + *H -> H2O
474  theReactionTable->SetReaction(1.44e10 * (1e-3 * m3 / (mole * s)), OH, H);
475  //------------------------------------------------------------------
476  // *H + *H -> H2
477  reactionData = new G4DNAMolecularReactionData(
478  1.20e10 * (1e-3 * m3 / (mole * s)), H, H);
479  reactionData->AddProduct(H2);
480  theReactionTable->SetReaction(reactionData);
481  //------------------------------------------------------------------
482  // H3O+ + OH- -> 2H2O
483  theReactionTable->SetReaction(1.43e11 * (1e-3 * m3 / (mole * s)), H3Op, OHm);
484  //------------------------------------------------------------------
485 }
486 
487 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
488 
490 {
492 
493  //===============================================================
494  // Extend vibrational to low energy
495  // Anyway, solvation of electrons is taken into account from 7.4 eV
496  // So below this threshold, for now, no accurate modeling is done
497  //
498  G4VProcess* process =
500  FindProcess("e-_G4DNAVibExcitation", "e-");
501 
502  if (process)
503  {
504  G4DNAVibExcitation* vibExcitation = (G4DNAVibExcitation*) process;
505  G4VEmModel* model = vibExcitation->EmModel();
506  G4DNASancheExcitationModel* sancheExcitationMod =
507  dynamic_cast<G4DNASancheExcitationModel*>(model);
508  if(sancheExcitationMod)
509  {
510  sancheExcitationMod->ExtendLowEnergyLimit(0.025 * eV);
511  }
512  }
513 
514  //===============================================================
515  // Modify elastic scattering models to avoid killing electrons
516  // at low energy
517  //
518  process =
520  FindProcess("e-_G4DNAElastic", "e-");
521 
522  if (process)
523  {
524  G4DNAElastic* vibExcitation = (G4DNAElastic*) process;
525  G4VEmModel* model = vibExcitation->EmModel();
526 
527  if(G4DNAChampionElasticModel* championMod =
528  dynamic_cast<G4DNAChampionElasticModel*>(model))
529  {
530  championMod->SetKillBelowThreshold(-1);
531  }
532  else if(G4DNAScreenedRutherfordElasticModel* screenRutherfordMod =
533  dynamic_cast<G4DNAScreenedRutherfordElasticModel*>(model))
534  {
535  screenRutherfordMod->SetKillBelowThreshold(-1);
536  }
537  else if (G4DNAUeharaScreenedRutherfordElasticModel* ueharaScreenRutherfordMod =
538  dynamic_cast<G4DNAUeharaScreenedRutherfordElasticModel*>(model))
539  {
540  ueharaScreenRutherfordMod->SetKillBelowThreshold(-1);
541  }
542  }
543 
544  //===============================================================
545  // *** Electron Solvatation ***
546  //
547  ph->RegisterProcess(
548  new G4DNAElectronSolvatation("e-_G4DNAElectronSolvatation"),
550 
551  //===============================================================
552  // Define processes for molecules
553  //
554  G4MoleculeTable* theMoleculeTable = G4MoleculeTable::Instance();
556  theMoleculeTable->GetDefintionIterator();
557  iterator.reset();
558  while (iterator())
559  {
560  G4MoleculeDefinition* moleculeDef = iterator.value();
561 
562  if (moleculeDef != G4H2O::Definition())
563  {
564  // G4cout << "Brownian motion added for : "
565  // << moleculeDef->GetName() << G4endl;
567  // brown->SetVerboseLevel(4);
568  ph->RegisterProcess(brown, moleculeDef);
569  }
570  else
571  {
572  moleculeDef->GetProcessManager()
574  G4DNAMolecularDissociation* dissociationProcess =
575  new G4DNAMolecularDissociation("H2O_DNAMolecularDecay");
576  dissociationProcess->SetDecayDisplacer(
577  moleculeDef, new G4DNAWaterDissociationDisplacer);
578  dissociationProcess->SetVerboseLevel(1);
579 // ph->RegisterProcess(dissociationProcess, moleculeDef);
580 
581  moleculeDef->GetProcessManager()
582  ->AddRestProcess(dissociationProcess, 1);
583  }
584  /*
585  * Warning : end of particles and processes are needed by
586  * EM Physics builders
587  */
588  }
589 
591 }
592 
593 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
594 
596  reactionTable)
597 {
598 
599  //=========================================
600  // Diffusion controlled reaction model
601  //=========================================
607  G4VDNAReactionModel* reactionRadiusComputer =
609  reactionTable->PrintTable(reactionRadiusComputer);
610 
616  G4DNAMolecularStepByStepModel* stepByStep =
618  stepByStep->SetReactionModel(reactionRadiusComputer);
619 // ((G4DNAMoleculeEncounterStepper*) stepByStep->GetTimeStepper())->
620 // SetVerbose(5);
621 
622  RegisterTimeStepModel(stepByStep, 0);
623 }
virtual void ConstructDissociationChannels()
static G4DLLIMPORT const DisplacementType AutoIonisation
static G4DLLIMPORT const DisplacementType B1A1_DissociationDecay
static G4Electron_aq * Definition()
static G4DLLIMPORT const DisplacementType NoDisplacement
float c_squared
Definition: hepunit.py:258
static G4H2O * Definition()
Definition: G4H2O.cc:46
virtual void ConstructReactionTable(G4DNAMolecularReactionTable *reactionTable)
static G4DLLIMPORT const DisplacementType A1B1_DissociationDecay
static G4DLLIMPORT const DisplacementType Ionisation_DissociationDecay
static const double m3
Definition: G4SIunits.hh:130
void AddDecayChannel(const G4MolecularConfiguration *molConf, const G4MolecularDissociationChannel *channel)
virtual void ConstructProcess()
virtual void ConstructTimeStepModel(G4DNAMolecularReactionTable *reactionTable)
static G4Electron * Definition()
Definition: G4Electron.cc:49
virtual ~G4EmDNAChemistry()
G4ProcessManager * GetProcessManager() const
float Avogadro
Definition: hepunit.py:253
G4VEmModel * EmModel(G4int index=1) const
static const double s
Definition: G4SIunits.hh:168
void SetReaction(G4double observedReactionRate, G4MolecularConfiguration *reactive1, G4MolecularConfiguration *reactive2)
const G4ElectronOccupancy * GetGroundStateElectronOccupancy() const
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5206
void SetChemistryList(G4VUserChemistryList *)
virtual void ConstructMolecule()
void AddProduct(const G4Molecule *, G4double=0)
static const double m2
Definition: G4SIunits.hh:129
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmDNAChemistry)
static G4H3O * Definition()
Definition: G4H3O.cc:47
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
G4int AddElectron(G4int orbit, G4int number=1)
static G4MoleculeTable * Instance()
void PrintTable(G4VDNAReactionModel *=0)
void AddProduct(G4MolecularConfiguration *molecule)
G4MolecularConfiguration * CreateConfiguration(const G4String &userIdentifier, const G4MoleculeDefinition *molDef, const G4String &configurationLabel, const G4ElectronOccupancy &eOcc)
static G4DNAChemistryManager * Instance()
static const double eV
Definition: G4SIunits.hh:212
G4MoleculeDefinitionIterator GetDefintionIterator()
static G4PhysicsListHelper * GetPhysicsListHelper()
static const double mole
Definition: G4SIunits.hh:283
static G4H2O2 * Definition()
Definition: G4H2O2.cc:46
static G4OH * Definition()
Definition: G4OH.cc:46
G4int AddRestProcess(G4VProcess *aProcess, G4int ord=ordDefault)
void SetDecayDisplacer(const G4ParticleDefinition *, G4VMolecularDecayDisplacer *)
G4MolecularConfiguration * NewConfigurationWithElectronOccupancy(const G4String &excitedStateLabel, const G4ElectronOccupancy &, double decayTime=0.)
static G4H2 * Definition()
Definition: G4H2.cc:46
void SetReactionModel(G4VDNAReactionModel *)
static G4ProcessTable * GetProcessTable()
G4int RemoveElectron(G4int orbit, G4int number=1)
void RegisterTimeStepModel(G4VITStepModel *timeStepModel, double startingTime=0)
static G4Hydrogen * Definition()
Definition: G4Hydrogen.cc:46
G4MolecularConfiguration * GetConfiguration(const G4String &, bool mustExist=true)