Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes
G4RadioactivationMessenger Class Reference

#include <G4RadioactivationMessenger.hh>

Inheritance diagram for G4RadioactivationMessenger:
G4UImessenger

Public Member Functions

G4bool CommandsShouldBeInMaster () const
 
 G4RadioactivationMessenger (G4Radioactivation *theRadioactivationContainer)
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
G4bool operator!= (const G4UImessenger &messenger) const
 
G4bool operator== (const G4UImessenger &messenger) const
 
void SetNewValue (G4UIcommand *command, G4String newValues)
 
 ~G4RadioactivationMessenger ()
 

Protected Member Functions

void AddUIcommand (G4UIcommand *newCommand)
 
G4String BtoS (G4bool b)
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
G4String DtoS (G4double a)
 
G4String ItoS (G4int i)
 
G4bool StoB (G4String s)
 
G4double StoD (G4String s)
 
G4int StoI (G4String s)
 
G4long StoL (G4String s)
 

Protected Attributes

G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Private Attributes

G4UIcmdWithABoolanaloguemcCmd
 
G4UIcmdWithABoolbrbiasCmd
 
G4UIcmdWithAStringdecaybiasprofileCmd
 
G4UIcmdWithADoubleAndUnithlthCmd
 
G4UIdirectoryrdmDirectory
 
G4UIcmdWithAStringsourcetimeprofileCmd
 
G4UIcmdWithAnIntegersplitnucleiCmd
 
G4RadioactivationtheRadioactivationContainer
 

Detailed Description

Definition at line 53 of file G4RadioactivationMessenger.hh.

Constructor & Destructor Documentation

◆ G4RadioactivationMessenger()

G4RadioactivationMessenger::G4RadioactivationMessenger ( G4Radioactivation theRadioactivationContainer)

Definition at line 42 of file G4RadioactivationMessenger.cc.

43 :theRadioactivationContainer(theRadioactivationContainer1)
44{
45 rdmDirectory = new G4UIdirectory("/process/had/rdm/");
46 rdmDirectory->SetGuidance("Controls the biased version of radioactive decay");
47
48 // Command to turn on/off variance reduction options
49 analoguemcCmd = new G4UIcmdWithABool("/process/had/rdm/analogueMC",this);
50 analoguemcCmd->SetGuidance("false: variance reduction method; true: analogue method");
51 analoguemcCmd->SetParameterName("AnalogueMC",true);
53
54 // Command to use branching ratio biasing or not
55 brbiasCmd = new G4UIcmdWithABool("/process/had/rdm/BRbias",this);
56 brbiasCmd->SetGuidance("false: no biasing; true: all branches are treated as equal");
57 brbiasCmd->SetParameterName("BRBias",true);
59
60 // Command to set the half-life thresold for isomer production
61 hlthCmd = new G4UIcmdWithADoubleAndUnit("/process/had/rdm/hlThreshold",this);
62 hlthCmd->SetGuidance("Set the h-l threshold for isomer production");
63 hlthCmd->SetParameterName("hlThreshold",false);
64 hlthCmd->SetUnitCategory("Time");
65
66 // Command to define the incident particle source time profile
67 sourcetimeprofileCmd = new G4UIcmdWithAString("/process/had/rdm/sourceTimeProfile",this);
69 ("Supply the name of the ascii file containing the source particle time profile");
70 sourcetimeprofileCmd->SetParameterName("STimeProfile",true);
72
73 // Command to define the incident particle source time profile
74 decaybiasprofileCmd = new G4UIcmdWithAString("/process/had/rdm/decayBiasProfile",this);
76 ("Supply the name of the ascii file containing the decay bias time profile");
77 decaybiasprofileCmd->SetParameterName("DBiasProfile",true);
79
80 // Command to set nuclei splitting parameter
81 splitnucleiCmd = new G4UIcmdWithAnInteger("/process/had/rdm/splitNuclei",this);
82 splitnucleiCmd->SetGuidance("Set number of splitting for the isotopes.");
83 splitnucleiCmd->SetParameterName("NSplit",true);
85 splitnucleiCmd->SetRange("NSplit>=1");
86}
G4Radioactivation * theRadioactivationContainer
G4UIcmdWithADoubleAndUnit * hlthCmd
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(G4bool defVal)
void SetUnitCategory(const char *unitCategory)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(const char *defVal)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(G4int defVal)
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
void SetRange(const char *rs)
Definition: G4UIcommand.hh:120

References analoguemcCmd, brbiasCmd, decaybiasprofileCmd, hlthCmd, rdmDirectory, G4UIcmdWithAString::SetDefaultValue(), G4UIcmdWithABool::SetDefaultValue(), G4UIcmdWithAnInteger::SetDefaultValue(), G4UIcommand::SetGuidance(), G4UIcmdWithABool::SetParameterName(), G4UIcmdWithADoubleAndUnit::SetParameterName(), G4UIcmdWithAnInteger::SetParameterName(), G4UIcmdWithAString::SetParameterName(), G4UIcommand::SetRange(), G4UIcmdWithADoubleAndUnit::SetUnitCategory(), sourcetimeprofileCmd, and splitnucleiCmd.

◆ ~G4RadioactivationMessenger()

G4RadioactivationMessenger::~G4RadioactivationMessenger ( )

Definition at line 89 of file G4RadioactivationMessenger.cc.

90{
91 delete rdmDirectory;
92 delete analoguemcCmd;
95 delete brbiasCmd;
96 delete splitnucleiCmd;
97 delete hlthCmd;
98}

References analoguemcCmd, brbiasCmd, decaybiasprofileCmd, hlthCmd, rdmDirectory, sourcetimeprofileCmd, and splitnucleiCmd.

Member Function Documentation

◆ AddUIcommand()

void G4UImessenger::AddUIcommand ( G4UIcommand newCommand)
protectedinherited

Definition at line 149 of file G4UImessenger.cc.

150{
151 G4cerr << "Warning : Old style definition of G4UIcommand <"
152 << newCommand->GetCommandPath() << ">." << G4endl;
153}
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
const G4String & GetCommandPath() const
Definition: G4UIcommand.hh:136

References G4cerr, G4endl, and G4UIcommand::GetCommandPath().

◆ BtoS()

G4String G4UImessenger::BtoS ( G4bool  b)
protectedinherited

Definition at line 98 of file G4UImessenger.cc.

99{
100 G4String vl = "0";
101 if(b)
102 vl = "true";
103 return vl;
104}

◆ CommandsShouldBeInMaster()

G4bool G4UImessenger::CommandsShouldBeInMaster ( ) const
inlineinherited

Definition at line 77 of file G4UImessenger.hh.

78 {
80 }
G4bool commandsShouldBeInMaster

References G4UImessenger::commandsShouldBeInMaster.

Referenced by G4UIcommand::G4UIcommandCommonConstructorCode().

◆ CreateCommand()

template<typename T >
T * G4UImessenger::CreateCommand ( const G4String cname,
const G4String dsc 
)
protectedinherited

Definition at line 110 of file G4UImessenger.hh.

111{
112 G4String path;
113 if(cname[0] != '/')
114 {
115 path = baseDirName + cname;
116 if(path[0] != '/')
117 path = "/" + path;
118 }
119
120 T* command = new T(path.c_str(), this);
121 command->SetGuidance(dsc.c_str());
122
123 return command;
124}
G4String baseDirName

References G4UImessenger::baseDirName.

◆ CreateDirectory()

void G4UImessenger::CreateDirectory ( const G4String path,
const G4String dsc,
G4bool  commandsToBeBroadcasted = true 
)
protectedinherited

Definition at line 156 of file G4UImessenger.cc.

158{
160
161 G4String fullpath = path;
162 if(fullpath.back() != '/')
163 fullpath.append("/");
164
165 G4UIcommandTree* tree = ui->GetTree()->FindCommandTree(fullpath.c_str());
166 if(tree != nullptr)
167 {
168 baseDirName = tree->GetPathName();
169 }
170 else
171 {
172 baseDir = new G4UIdirectory(fullpath.c_str(), commandsToBeBroadcasted);
173 baseDirName = fullpath;
174 baseDir->SetGuidance(dsc.c_str());
175 }
176}
const G4String & GetPathName() const
G4UIcommandTree * FindCommandTree(const char *commandPath)
G4UIcommandTree * GetTree() const
Definition: G4UImanager.hh:186
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
G4UIdirectory * baseDir

References G4UImessenger::baseDir, G4UImessenger::baseDirName, G4UIcommandTree::FindCommandTree(), G4UIcommandTree::GetPathName(), G4UImanager::GetTree(), G4UImanager::GetUIpointer(), and G4UIcommand::SetGuidance().

Referenced by G4MoleculeShootMessenger::G4MoleculeShootMessenger(), and G4UImessenger::G4UImessenger().

◆ DtoS()

G4String G4UImessenger::DtoS ( G4double  a)
protectedinherited

Definition at line 90 of file G4UImessenger.cc.

91{
92 std::ostringstream os;
93 os << a;
94 return G4String(os.str());
95}

Referenced by G4ScoreQuantityMessenger::FilterCommands(), and G4UIcontrolMessenger::SetNewValue().

◆ GetCurrentValue()

G4String G4UImessenger::GetCurrentValue ( G4UIcommand command)
virtualinherited

Reimplemented in G4ScoreQuantityMessenger, G4VisCommandModelCreate< Factory >, G4VisCommandListManagerList< Manager >, G4VisCommandListManagerSelect< Manager >, G4VisCommandManagerMode< Manager >, G4ToolsAnalysisMessenger, G4ScoringMessenger, G4EvManMessenger, G4GeneralParticleSourceMessenger, G4ParticleGunMessenger, G4GeometryMessenger, G4GenericMessenger, G4UIcontrolMessenger, GFlashShowerModelMessenger, G4DecayTableMessenger, G4ParticleMessenger, G4ParticlePropertyMessenger, G4tgrMessenger, G4PersistencyCenterMessenger, G4ProductionCutsTableMessenger, G4SchedulerMessenger, G4VITSteppingVerbose, G4MoleculeShootMessenger, G4MoleculeGunMessenger, G4ProcessManagerMessenger, G4ProcessTableMessenger, G4MatScanMessenger, G4RunMessenger, G4UserPhysicsListMessenger, G4TrackingMessenger, G4GMocrenMessenger, G4HepRepMessenger, G4VisCommandAbortReviewKeptEvents, G4VisCommandDrawOnlyToBeKeptEvents, G4VisCommandEnable, G4VisCommandList, G4VisCommandReviewKeptEvents, G4VisCommandVerbose, G4VisCommandGeometryList, G4VisCommandGeometryRestore, G4VisCommandGeometrySetColour, G4VisCommandGeometrySetDaughtersInvisible, G4VisCommandGeometrySetForceAuxEdgeVisible, G4VisCommandGeometrySetForceCloud, G4VisCommandGeometrySetForceSolid, G4VisCommandGeometrySetForceLineSegmentsPerCircle, G4VisCommandGeometrySetForceWireframe, G4VisCommandGeometrySetLineStyle, G4VisCommandGeometrySetLineWidth, G4VisCommandGeometrySetVisibility, G4VisCommandSceneActivateModel, G4VisCommandSceneCreate, G4VisCommandSceneEndOfEventAction, G4VisCommandSceneEndOfRunAction, G4VisCommandSceneList, G4VisCommandSceneNotifyHandlers, G4VisCommandSceneRemoveModel, G4VisCommandSceneSelect, G4VisCommandSceneShowExtents, G4VisCommandSceneAddArrow, G4VisCommandSceneAddArrow2D, G4VisCommandSceneAddAxes, G4VisCommandSceneAddDate, G4VisCommandSceneAddDigis, G4VisCommandSceneAddEventID, G4VisCommandSceneAddExtent, G4VisCommandSceneAddElectricField, G4VisCommandSceneAddFrame, G4VisCommandSceneAddGPS, G4VisCommandSceneAddGhosts, G4VisCommandSceneAddHits, G4VisCommandSceneAddLine, G4VisCommandSceneAddLine2D, G4VisCommandSceneAddLocalAxes, G4VisCommandSceneAddLogicalVolume, G4VisCommandSceneAddLogo, G4VisCommandSceneAddLogo2D, G4VisCommandSceneAddMagneticField, G4VisCommandSceneAddPSHits, G4VisCommandSceneAddScale, G4VisCommandSceneAddText, G4VisCommandSceneAddText2D, G4VisCommandSceneAddTrajectories, G4VisCommandSceneAddUserAction, G4VisCommandSceneAddVolume, G4VisCommandSceneAddPlotter, G4VisCommandSceneHandlerAttach, G4VisCommandSceneHandlerCreate, G4VisCommandSceneHandlerList, G4VisCommandSceneHandlerSelect, G4VisCommandSetArrow3DLineSegmentsPerCircle, G4VisCommandSetColour, G4VisCommandSetExtentForField, G4VisCommandSetLineWidth, G4VisCommandSetTextColour, G4VisCommandSetTextLayout, G4VisCommandSetTextSize, G4VisCommandSetTouchable, G4VisCommandSetVolumeForField, G4VisCommandsTouchable, G4VisCommandsTouchableSet, G4VisCommandViewerAddCutawayPlane, G4VisCommandViewerCentreOn, G4VisCommandViewerChangeCutawayPlane, G4VisCommandViewerClear, G4VisCommandViewerClearCutawayPlanes, G4VisCommandViewerClearTransients, G4VisCommandViewerClearVisAttributesModifiers, G4VisCommandViewerClone, G4VisCommandViewerColourByDensity, G4VisCommandViewerCopyViewFrom, G4VisCommandViewerCreate, G4VisCommandViewerDolly, G4VisCommandViewerFlush, G4VisCommandViewerInterpolate, G4VisCommandViewerList, G4VisCommandViewerPan, G4VisCommandViewerReset, G4VisCommandViewerRefresh, G4VisCommandViewerRebuild, G4VisCommandViewerSave, G4VisCommandViewerScale, G4VisCommandViewerSelect, G4VisCommandViewerUpdate, G4VisCommandViewerZoom, G4VisCommandViewerDefaultHiddenEdge, G4VisCommandViewerDefaultStyle, G4VisCommandsViewerSet, G4VModelCommand< T >, G4VModelCommand< M >, G4RTMessenger, G4ASCIITreeMessenger, G4VtkMessenger, G4PolarizationMessenger, and G4DNAChemistryManager.

Definition at line 58 of file G4UImessenger.cc.

59{
60 G4String nullString;
61 return nullString;
62}

Referenced by G4UIcommand::DoIt(), and G4UIcommand::GetCurrentValue().

◆ ItoS()

G4String G4UImessenger::ItoS ( G4int  i)
protectedinherited

Definition at line 82 of file G4UImessenger.cc.

83{
84 std::ostringstream os;
85 os << i;
86 return G4String(os.str());
87}

Referenced by G4GenericMessenger::DeclareMethod(), and G4ParticleGunMessenger::GetCurrentValue().

◆ operator!=()

G4bool G4UImessenger::operator!= ( const G4UImessenger messenger) const
inherited

Definition at line 76 of file G4UImessenger.cc.

77{
78 return this != &messenger;
79}

◆ operator==()

G4bool G4UImessenger::operator== ( const G4UImessenger messenger) const
inherited

Definition at line 70 of file G4UImessenger.cc.

71{
72 return this == &messenger;
73}

◆ SetNewValue()

void G4RadioactivationMessenger::SetNewValue ( G4UIcommand command,
G4String  newValues 
)
virtual

Reimplemented from G4UImessenger.

Definition at line 101 of file G4RadioactivationMessenger.cc.

102{
103 if ( command == analoguemcCmd ) { theRadioactivationContainer->
104 SetAnalogueMonteCarlo( analoguemcCmd->GetNewBoolValue( newValues ) );
105 } else if ( command == brbiasCmd ) { theRadioactivationContainer->
106 SetBRBias( brbiasCmd->GetNewBoolValue( newValues ) );
107 } else if ( command == sourcetimeprofileCmd ) { theRadioactivationContainer->
108 SetSourceTimeProfile( newValues );
109 } else if ( command == decaybiasprofileCmd ) { theRadioactivationContainer->
110 SetDecayBias( newValues );
111 } else if ( command == splitnucleiCmd ) { theRadioactivationContainer->
112 SetSplitNuclei( splitnucleiCmd->GetNewIntValue( newValues ) );
113 } else if ( command == hlthCmd ) { theRadioactivationContainer->
114 SetHLThreshold( hlthCmd->GetNewDoubleValue( newValues ) );
115 }
116}
static G4bool GetNewBoolValue(const char *paramString)
static G4double GetNewDoubleValue(const char *paramString)
static G4int GetNewIntValue(const char *paramString)

References analoguemcCmd, brbiasCmd, decaybiasprofileCmd, G4UIcmdWithABool::GetNewBoolValue(), G4UIcmdWithADoubleAndUnit::GetNewDoubleValue(), G4UIcmdWithAnInteger::GetNewIntValue(), hlthCmd, sourcetimeprofileCmd, splitnucleiCmd, and theRadioactivationContainer.

◆ StoB()

G4bool G4UImessenger::StoB ( G4String  s)
protectedinherited

Definition at line 137 of file G4UImessenger.cc.

138{
140 G4bool vl = false;
141 if(v == "Y" || v == "YES" || v == "1" || v == "T" || v == "TRUE")
142 {
143 vl = true;
144 }
145 return vl;
146}
bool G4bool
Definition: G4Types.hh:86
G4String to_upper_copy(G4String str)
Return uppercase copy of string.

References G4StrUtil::to_upper_copy().

Referenced by G4LocalThreadCoutMessenger::SetNewValue(), G4CascadeParamMessenger::SetNewValue(), G4ScoreQuantityMessenger::SetNewValue(), and G4ScoringMessenger::SetNewValue().

◆ StoD()

G4double G4UImessenger::StoD ( G4String  s)
protectedinherited

◆ StoI()

G4int G4UImessenger::StoI ( G4String  s)
protectedinherited

◆ StoL()

G4long G4UImessenger::StoL ( G4String  s)
protectedinherited

Definition at line 117 of file G4UImessenger.cc.

118{
119 G4long vl;
120 const char* t = str;
121 std::istringstream is(t);
122 is >> vl;
123 return vl;
124}
long G4long
Definition: G4Types.hh:87

Referenced by G4RunMessenger::SetNewValue().

Field Documentation

◆ analoguemcCmd

G4UIcmdWithABool* G4RadioactivationMessenger::analoguemcCmd
private

◆ baseDir

G4UIdirectory* G4UImessenger::baseDir = nullptr
protectedinherited

◆ baseDirName

G4String G4UImessenger::baseDirName = ""
protectedinherited

◆ brbiasCmd

G4UIcmdWithABool* G4RadioactivationMessenger::brbiasCmd
private

◆ commandsShouldBeInMaster

G4bool G4UImessenger::commandsShouldBeInMaster = false
protectedinherited

◆ decaybiasprofileCmd

G4UIcmdWithAString* G4RadioactivationMessenger::decaybiasprofileCmd
private

◆ hlthCmd

G4UIcmdWithADoubleAndUnit* G4RadioactivationMessenger::hlthCmd
private

◆ rdmDirectory

G4UIdirectory* G4RadioactivationMessenger::rdmDirectory
private

◆ sourcetimeprofileCmd

G4UIcmdWithAString* G4RadioactivationMessenger::sourcetimeprofileCmd
private

◆ splitnucleiCmd

G4UIcmdWithAnInteger* G4RadioactivationMessenger::splitnucleiCmd
private

◆ theRadioactivationContainer

G4Radioactivation* G4RadioactivationMessenger::theRadioactivationContainer
private

Definition at line 62 of file G4RadioactivationMessenger.hh.

Referenced by SetNewValue().


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