Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/radioactive_decay/src/G4RadioactivationMessenger.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 ////////////////////////////////////////////////////////////////////////////////
 27 //
 28 //  GEANT4 Class source file
 29 //
 30 //  G4RadioactivationMessenger
 31 //
 32 //  Author: D.H. Wright (SLAC)
 33 //  Date:   29 August 2017
 34 //
 35 //  Based on the code of F. Lei and P.R. Truscott.
 36 //
 37 ////////////////////////////////////////////////////////////////////////////////
 38 
 39 #include "G4RadioactivationMessenger.hh"
 40 #include "G4NuclearLevelData.hh"
 41 #include "G4RadioactiveDecay.hh"
 42 
 43 #include <sstream>
 44 #include "G4HadronicException.hh"
 45 
 46 
 47 G4RadioactivationMessenger::G4RadioactivationMessenger(G4RadioactiveDecay* ptr)
 48  :theRadDecay(ptr)
 49 {
 50   rdmDirectory = new G4UIdirectory("/process/had/rdm/");
 51   rdmDirectory->SetGuidance("Controls the biased version of radioactive decay");
 52 
 53   // Command to turn on/off variance reduction options
 54   analoguemcCmd = new G4UIcmdWithABool("/process/had/rdm/analogueMC",this);
 55   analoguemcCmd->SetGuidance("false: variance reduction method; true: analogue method");
 56   analoguemcCmd->SetParameterName("AnalogueMC",true);
 57   analoguemcCmd->SetDefaultValue(true);
 58   
 59   // Command to use branching ratio biasing or not
 60   brbiasCmd = new G4UIcmdWithABool("/process/had/rdm/BRbias",this);
 61   brbiasCmd->SetGuidance("false: no biasing; true: all branches are treated as equal");
 62   brbiasCmd->SetParameterName("BRBias",true);
 63   brbiasCmd->SetDefaultValue(true);
 64   
 65   // Command to set the half-life thresold for isomer production
 66   hlthCmd = new G4UIcmdWithADoubleAndUnit("/process/had/rdm/hlThreshold",this);
 67   hlthCmd->SetGuidance("Set the h-l threshold for isomer production");
 68   hlthCmd->SetParameterName("hlThreshold",false);
 69   hlthCmd->SetUnitCategory("Time");
 70   
 71   // Command to define the incident particle source time profile
 72   sourcetimeprofileCmd = new G4UIcmdWithAString("/process/had/rdm/sourceTimeProfile",this);
 73   sourcetimeprofileCmd->SetGuidance 
 74     ("Supply the name of the ascii file containing the source particle time profile");
 75   sourcetimeprofileCmd->SetParameterName("STimeProfile",true);
 76   sourcetimeprofileCmd->SetDefaultValue("source.data");
 77   
 78   // Command to define the incident particle source time profile
 79   decaybiasprofileCmd = new G4UIcmdWithAString("/process/had/rdm/decayBiasProfile",this);
 80   decaybiasprofileCmd->SetGuidance 
 81     ("Supply the name of the ascii file containing the decay bias time profile");
 82   decaybiasprofileCmd->SetParameterName("DBiasProfile",true);
 83   decaybiasprofileCmd->SetDefaultValue("bias.data");
 84   
 85   // Command to set nuclei splitting parameter
 86   splitnucleiCmd = new G4UIcmdWithAnInteger("/process/had/rdm/splitNuclei",this);
 87   splitnucleiCmd->SetGuidance("Set number of splitting for the isotopes.");
 88   splitnucleiCmd->SetParameterName("NSplit",true);
 89   splitnucleiCmd->SetDefaultValue(1);
 90   splitnucleiCmd->SetRange("NSplit>=1");
 91 }
 92 
 93 
 94 G4RadioactivationMessenger::~G4RadioactivationMessenger()
 95 {
 96   delete rdmDirectory;
 97   delete analoguemcCmd;
 98   delete sourcetimeprofileCmd;
 99   delete decaybiasprofileCmd;
100   delete brbiasCmd;
101   delete splitnucleiCmd;
102   delete hlthCmd;
103 }
104 
105 
106 void G4RadioactivationMessenger::SetNewValue(G4UIcommand* command, G4String newValues)
107 {
108   if ( command == analoguemcCmd ) { theRadDecay->
109     SetAnalogueMonteCarlo( analoguemcCmd->GetNewBoolValue( newValues ) );
110   } else if ( command == brbiasCmd ) { theRadDecay->
111     SetBRBias( brbiasCmd->GetNewBoolValue( newValues ) );
112   } else if ( command == sourcetimeprofileCmd ) { theRadDecay->
113     SetSourceTimeProfile( newValues );
114   } else if ( command == decaybiasprofileCmd ) { theRadDecay->
115     SetDecayBias( newValues );
116   } else if ( command == splitnucleiCmd ) { theRadDecay->
117     SetSplitNuclei( splitnucleiCmd->GetNewIntValue( newValues ) );
118   } else if ( command == hlthCmd ) { theRadDecay->
119     SetHLThreshold( hlthCmd->GetNewDoubleValue( newValues ) );
120   }
121 }
122 
123