Geant4-11
G4GlobalMagFieldMessenger.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// class G4GlobalMagFieldMessenger
27//
28// Implementation of the G4GlobalMagFieldMessenger class
29//
30// Author: Ivana Hrivnacova, 28/08/2013 (ivana@ipno.in2p3.fr)
31// --------------------------------------------------------------------
32
34#include "G4UniformMagField.hh"
35#include "G4FieldManager.hh"
37#include "G4UIdirectory.hh"
40#include "G4UnitsTable.hh"
41#include "G4SystemOfUnits.hh"
42
43//______________________________________________________________________________
44
47{
48 fDirectory = new G4UIdirectory("/globalField/");
49 fDirectory->SetGuidance("Global uniform magnetic field UI commands");
50
51 fSetValueCmd = new G4UIcmdWith3VectorAndUnit("/globalField/setValue",this);
52 fSetValueCmd->SetGuidance("Set uniform magnetic field value.");
53 fSetValueCmd->SetParameterName("Bx", "By", "By", false);
54 fSetValueCmd->SetUnitCategory("Magnetic flux density");
56
57 fSetVerboseCmd = new G4UIcmdWithAnInteger("/globalField/verbose",this);
58 fSetVerboseCmd->SetGuidance("Set verbose level: ");
59 fSetVerboseCmd->SetGuidance(" 0: no output");
60 fSetVerboseCmd->SetGuidance(" 1: printing new field value");
61 fSetVerboseCmd->SetParameterName("globalFieldVerbose", false);
62 fSetVerboseCmd->SetRange("globalFieldVerbose>=0");
64
65 // Create field
66 fMagField = new G4UniformMagField(value);
67
68 // Set field value (the field is not activated if value is zero)
69 SetField(value, "G4GlobalMagFieldMessenger::G4GlobalMagFieldMessenger");
70}
71
72//______________________________________________________________________________
73
75{
76 delete fMagField;
77 delete fSetValueCmd;
78 delete fSetVerboseCmd;
79 delete fDirectory;
80}
81
82//______________________________________________________________________________
83
85 const G4String& /*inFunction*/)
86{
87 // Get field manager (or create it if it does not yet exist)
88 G4FieldManager* fieldManager
90
91 // Inactivate field if its value is zero
92 if ( value == G4ThreeVector() )
93 {
94 fieldManager->SetDetectorField(0);
95 fieldManager->CreateChordFinder(0);
96
97 if ( fVerboseLevel > 0 )
98 {
99 G4cout << "Magnetic field is inactive, fieldValue = (0,0,0)." << G4endl;
100 }
101 }
102 else
103 {
104 fMagField->SetFieldValue(value);
105 fieldManager->SetDetectorField(fMagField);
106 fieldManager->CreateChordFinder(fMagField);
107
108 if ( fVerboseLevel > 0 )
109 {
110 G4cout << "Magnetic field is active, fieldValue = ("
111 << G4BestUnit(value, "Magnetic flux density")
112 << ")." << G4endl;
113 }
114 }
115}
116
117//______________________________________________________________________________
118
120{
121 if ( command == fSetValueCmd )
122 {
124 "G4GlobalMagFieldMessenger::SetNewValue");
125 }
126 else if ( command == fSetVerboseCmd )
127 {
129 }
130}
131
132//______________________________________________________________________________
133
135{
136 SetField(value, "G4GlobalMagFieldMessenger::SetFieldValue");
137}
138
139//______________________________________________________________________________
140
142{
144
145 return G4ThreeVector();
146}
@ G4State_Idle
@ G4State_PreInit
#define G4BestUnit(a, b)
CLHEP::Hep3Vector G4ThreeVector
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4bool SetDetectorField(G4Field *detectorField, G4int failMode=0)
void CreateChordFinder(G4MagneticField *detectorMagField)
void SetVerboseLevel(G4int verboseLevel)
void SetField(const G4ThreeVector &value, const G4String &inFunction)
G4GlobalMagFieldMessenger(const G4ThreeVector &value=G4ThreeVector())
G4UIcmdWith3VectorAndUnit * fSetValueCmd
G4UIcmdWithAnInteger * fSetVerboseCmd
void SetFieldValue(const G4ThreeVector &value)
virtual void SetNewValue(G4UIcommand *, G4String)
static G4TransportationManager * GetTransportationManager()
G4FieldManager * GetFieldManager() const
static G4ThreeVector GetNew3VectorValue(const char *paramString)
void SetParameterName(const char *theNameX, const char *theNameY, const char *theNameZ, G4bool omittable, G4bool currentAsDefault=false)
void SetUnitCategory(const char *unitCategory)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
static G4int GetNewIntValue(const char *paramString)
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
void SetRange(const char *rs)
Definition: G4UIcommand.hh:120
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:288
G4ThreeVector GetConstantFieldValue() const
void SetFieldValue(const G4ThreeVector &newFieldValue)