Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4BlineTracerMessenger.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 /// \file field/BlineTracer/src/G4BlineTracerMessenger.cc
27 /// \brief Implementation of the G4BlineTracerMessenger class
28 //
29 //
30 // $Id: G4BlineTracerMessenger.cc 68021 2013-03-13 13:36:07Z gcosmo $
31 //
32 // --------------------------------------------------------------------
33 //
34 // G4BlineTracerMessenger implementation
35 //
36 // --------------------------------------------------------------------
37 // Author: Laurent Desorgher (desorgher@phim.unibe.ch)
38 // Created - 2003-10-06
39 // --------------------------------------------------------------------
40 
42 #include "G4BlineTracer.hh"
43 #include "G4BlineEventAction.hh"
44 
45 #include "G4UImessenger.hh"
46 #include "G4UIdirectory.hh"
47 #include "G4UIcommand.hh"
48 #include "G4UIcmdWithADouble.hh"
50 #include "G4UIcmdWithAnInteger.hh"
52 #include "G4UIcmdWith3Vector.hh"
53 #include "G4UIcmdWithABool.hh"
54 
55 ///////////////////////////////////////////////////////////////////////////
56 
58 {
59  fTheBlineTool = aBlineTool;
60  fBlineToolDir = new G4UIdirectory("/vis/blineTracer/");
61  fBlineToolDir->SetGuidance("Commands to trace and visualise magnetic field lines.");
62  fBlineToolDir->SetGuidance("These commands work only if a magnetic-field is set");
63  fBlineToolDir->SetGuidance("in the application.");
64 
65  // commands
66 
67  fBlineCmd = new G4UIcmdWithAnInteger("/vis/blineTracer/computeBline",this);
68  fBlineCmd->SetGuidance("Compute magnetic field lines for visualisation.");
69  fBlineCmd->SetParameterName("nb_of_lines",false);
71 
72  fSetMaxTrackingStepCmd =
73  new G4UIcmdWithADoubleAndUnit("/vis/blineTracer/setMaxStepLength",this);
74  fSetMaxTrackingStepCmd->SetGuidance("Set the maximum length of tracking step");
75  fSetMaxTrackingStepCmd->SetGuidance("when integrating magnetic field line.");
76  fSetMaxTrackingStepCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
77 
78  fSetDrawColourCmd = new G4UIcmdWith3Vector("/vis/blineTracer/setColour",this);
79  fSetDrawColourCmd->SetGuidance("Set the colour drawing trajectories");
80  fSetDrawColourCmd->SetGuidance("and magnetic field lines.");
82 
83  fSetDrawBlineCmd = new G4UIcmdWithABool("/vis/blineTracer/stockLines",this);
84  fSetDrawBlineCmd->SetGuidance("If true field lines are stocked in lines");
85  fSetDrawBlineCmd->SetGuidance("to be drawn.");
86  fSetDrawBlineCmd->SetParameterName("StockLines",false);
88 
89  fSetDrawPointsCmd = new G4UIcmdWithABool("/vis/blineTracer/stockPoints",this);
90  fSetDrawPointsCmd->SetGuidance("If true step field line points are stocked");
91  fSetDrawPointsCmd->SetGuidance("in vector of points to be drawn.");
92  fSetDrawPointsCmd->SetParameterName("StockPoints",false);
94 
95  fSetPointSizeCmd = new G4UIcmdWithADouble("/vis/blineTracer/setPointSize",this);
96  fSetPointSizeCmd->SetGuidance("Set the size of points for drawing.");
97  fSetPointSizeCmd->SetParameterName("StepSize",false);
99 
100  fDrawCmd = new G4UIcmdWithoutParameter("/vis/blineTracer/show",this);
101  fDrawCmd->SetGuidance("Show the stored magnetic field lines.");
103 
104  fResetCmd =
105  new G4UIcmdWithoutParameter("/vis/blineTracer/resetMaterialToBeDrawn",this);
106  fResetCmd->SetGuidance("Clear the vectors of lines and points to be drawn.");
108 }
109 
110 ///////////////////////////////////////////////////////////////////////////
111 
113 {
114  delete fResetCmd;
115  delete fDrawCmd;
116  delete fSetPointSizeCmd;
117  delete fSetDrawPointsCmd;
118  delete fSetDrawBlineCmd;
119  delete fSetDrawColourCmd;
120  delete fSetMaxTrackingStepCmd;
121  delete fBlineCmd;
122  delete fBlineToolDir;
123 }
124 
125 ///////////////////////////////////////////////////////////////////////////
126 
128  G4String newValues )
129 {
130  if (command == fBlineCmd)
131  fTheBlineTool->ComputeBlines(1);
132  else if( command == fSetMaxTrackingStepCmd )
133  fTheBlineTool->SetMaxTrackingStep(fSetMaxTrackingStepCmd
134  ->GetNewDoubleValue(newValues));
135  else if( command == fSetDrawBlineCmd )
136  fTheBlineTool->GetEventAction()->SetDrawBline(fSetDrawBlineCmd
137  ->GetNewBoolValue(newValues));
138  else if( command == fSetDrawColourCmd )
139  {
140  G4ThreeVector vec=fSetDrawColourCmd->GetNew3VectorValue(newValues);
141  fTheBlineTool->GetEventAction()->
142  SetDrawColour(G4Colour(vec.x(),vec.y(),vec.z()));
143  }
144  else if( command == fSetDrawPointsCmd )
145  fTheBlineTool->GetEventAction()->SetDrawPoints(fSetDrawPointsCmd
146  ->GetNewBoolValue(newValues));
147  else if( command == fSetPointSizeCmd )
148  fTheBlineTool->GetEventAction()->SetPointSize(fSetPointSizeCmd
149  ->GetNewDoubleValue(newValues));
150  else if( command == fDrawCmd )
151  fTheBlineTool->GetEventAction()->DrawFieldLines(.5,45.,45.);
152  else if( command == fResetCmd )
153  fTheBlineTool->GetEventAction()->ResetVectorObjectToBeDrawn();
154 }
void SetDrawBline(G4bool aBool)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
double x() const
void DrawFieldLines(G4double zoom, G4double theta, G4double phi)
void SetDrawPoints(G4bool aBool)
double z() const
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
Definition of the G4BlineTracer class.
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:161
Definition of the G4BlineTracerMessenger class.
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:225
void ComputeBlines(G4int nlines)
void SetMaxTrackingStep(G4double max_step)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetPointSize(G4double aVal)
double y() const
G4BlineEventAction * GetEventAction()
Definition of the G4BlineEventAction class.
G4BlineTracerMessenger(G4BlineTracer *aBlineTool)
virtual void SetNewValue(G4UIcommand *command, G4String newValues)