Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VScoringMesh.hh
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 //
27 // $Id: G4VScoringMesh.hh 99154 2016-09-07 08:06:30Z gcosmo $
28 //
29 
30 #ifndef G4VScoringMesh_h
31 #define G4VScoringMesh_h 1
32 
33 #include "globals.hh"
34 #include "G4THitsMap.hh"
35 #include "G4RotationMatrix.hh"
36 #include "G4StatDouble.hh"
37 
38 class G4VPhysicalVolume;
39 class G4LogicalVolume;
41 class G4VPrimitiveScorer;
42 class G4VSDFilter;
43 class G4VScoreColorMap;
45 
46 #include <map>
47 
51 typedef std::map< G4String, RunScore* > MeshScoreMap;
52 // class description:
53 //
54 // This class represents a parallel world for interactive scoring purposes.
55 //
56 
58 {
59  public:
60  G4VScoringMesh(const G4String& wName);
61  virtual ~G4VScoringMesh();
62 
63  public: // with description
64  // a pure virtual function to construct this mesh geometry
65  void Construct(G4VPhysicalVolume* fWorldPhys);
66 
67  void WorkerConstruct(G4VPhysicalVolume* fWorldPhys);
68 
69  protected:
70  virtual void SetupGeometry(G4VPhysicalVolume * fWorldPhys) = 0;
71 
72  public: // with description
73  // list infomration of this mesh
74  virtual void List() const;
75 
76  public: // with description
77  // get the world name
78  inline const G4String& GetWorldName() const
79  { return fWorldName; }
80  // get whether this mesh is active or not
81  inline G4bool IsActive() const
82  { return fActive; }
83  // set an activity of this mesh
84  inline void Activate(G4bool vl = true)
85  { fActive = vl; }
86  // get the shape of this mesh
87  inline MeshShape GetShape() const
88  { return fShape; }
89  // accumulate hits in a registered primitive scorer
90  void Accumulate(G4THitsMap<G4double> * map);
92  // merge same kind of meshes
93  void Merge(const G4VScoringMesh * scMesh);
94  // dump information of primitive socrers registered in this mesh
95  void Dump();
96  // draw a projected quantity on a current viewer
97  void DrawMesh(const G4String& psName,G4VScoreColorMap* colorMap,G4int axflg=111);
98  // draw a column of a quantity on a current viewer
99  void DrawMesh(const G4String& psName,G4int idxPlane,G4int iColumn,G4VScoreColorMap* colorMap);
100  // draw a projected quantity on a current viewer
101  virtual void Draw(RunScore * map, G4VScoreColorMap* colorMap, G4int axflg=111) = 0;
102  // draw a column of a quantity on a current viewer
103  virtual void DrawColumn(RunScore * map, G4VScoreColorMap* colorMap,
104  G4int idxProj, G4int idxColumn) = 0;
105  // reset registered primitive scorers
106  void ResetScore();
107 
108  // set size of this mesh
109  void SetSize(G4double size[3]);
110  // get size of this mesh
111  G4ThreeVector GetSize() const;
112  // set position of center of this mesh
113  void SetCenterPosition(G4double centerPosition[3]);
114  // get position of center of this mesh
116  // set a rotation angle around the x axis
117  void RotateX(G4double delta);
118  // set a rotation angle around the y axis
119  void RotateY(G4double delta);
120  // set a rotation angle around the z axis
121  void RotateZ(G4double delta);
122  // get a rotation matrix
124  if(fRotationMatrix) return *fRotationMatrix;
125  else return G4RotationMatrix::IDENTITY;
126  }
127  // set number of segments of this mesh
128  void SetNumberOfSegments(G4int nSegment[3]);
129  // get number of segments of this mesh
130  void GetNumberOfSegments(G4int nSegment[3]);
131 
132  // register a primitive scorer to the MFD & set it to the current primitive scorer
134  // register a filter to a current primtive scorer
135  void SetFilter(G4VSDFilter * filter);
136  // set a primitive scorer to the current one by the name
138  // find registered primitive scorer by the name
139  G4bool FindPrimitiveScorer(const G4String & psname);
140  // get whether current primitive scorer is set or not
142  if(fCurrentPS == nullptr) return true;
143  else return false;
144  }
145  // get unit of primitive scorer by the name
146  G4String GetPSUnit(const G4String & psname);
147  // get unit of current primitive scorer
149  // set unit of current primitive scorer
150  void SetCurrentPSUnit(const G4String& unit);
151  // get unit value of primitive scorer by the name
152  G4double GetPSUnitValue(const G4String & psname);
153  // set PS name to be drawn
154  void SetDrawPSName(const G4String & psname) {fDrawPSName = psname;}
155 
156  // get axis names of the hierarchical division in the divided order
157  void GetDivisionAxisNames(G4String divisionAxisNames[3]);
158 
159  // set current primitive scorer to NULL
161  // set verbose level
162  inline void SetVerboseLevel(G4int vl)
163  { verboseLevel = vl; }
164  // get the primitive scorer map
165  inline MeshScoreMap GetScoreMap() const
166  { return fMap; }
167  // get whether this mesh setup has been ready
168  inline G4bool ReadyForQuantity() const
169  { return (sizeIsSet && nMeshIsSet); }
170 
171 protected:
172  // get registered primitive socrer by the name
174 
175 protected:
181 
186 
189 
191 
194 
198 
200 
202 
203 public:
205  { fMeshElementLogical = val; }
207  { return fMeshElementLogical; }
208 
209 protected:
212 public:
214  { fParallelWorldProcess = proc; }
216  { return fParallelWorldProcess; }
218  {
220  fMeshElementLogical = nullptr;
221  }
222 };
223 
224 #endif
225 
void GetNumberOfSegments(G4int nSegment[3])
const XML_Char * name
Definition: expat.h:151
const G4String & GetWorldName() const
G4VScoringMesh(const G4String &wName)
G4RotationMatrix GetRotationMatrix() const
void Activate(G4bool vl=true)
void DrawMesh(const G4String &psName, G4VScoreColorMap *colorMap, G4int axflg=111)
void SetParallelWorldProcess(G4ParallelWorldProcess *proc)
G4ThreeVector fCenterPosition
void Merge(const G4VScoringMesh *scMesh)
G4bool FindPrimitiveScorer(const G4String &psname)
void RotateY(G4double delta)
G4VPrimitiveScorer * GetPrimitiveScorer(const G4String &name)
void GeometryHasBeenDestroyed()
void SetPrimitiveScorer(G4VPrimitiveScorer *ps)
void SetMeshElementLogical(G4LogicalVolume *val)
void SetCurrentPrimitiveScorer(const G4String &name)
G4MultiFunctionalDetector * fMFD
void SetSize(G4double size[3])
G4double GetPSUnitValue(const G4String &psname)
void Accumulate(G4THitsMap< G4double > *map)
void RotateX(G4double delta)
G4THitsMap< G4StatDouble > RunScore
G4double fSize[3]
void SetFilter(G4VSDFilter *filter)
void WorkerConstruct(G4VPhysicalVolume *fWorldPhys)
G4ThreeVector GetSize() const
G4String GetPSUnit(const G4String &psname)
int G4int
Definition: G4Types.hh:78
void Construct(G4VPhysicalVolume *fWorldPhys)
G4ParallelWorldProcess * fParallelWorldProcess
void SetNullToCurrentPrimitiveScorer()
G4bool IsCurrentPrimitiveScorerNull()
static DLL_API const HepRotation IDENTITY
Definition: Rotation.h:369
virtual void Draw(RunScore *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0
G4double fDrawUnitValue
MeshScoreMap fMap
MeshShape GetShape() const
virtual void List() const
G4RotationMatrix * fRotationMatrix
bool G4bool
Definition: G4Types.hh:79
G4ThreeVector GetTranslation() const
G4THitsMap< G4double > EventScore
G4bool ReadyForQuantity() const
virtual void SetupGeometry(G4VPhysicalVolume *fWorldPhys)=0
void SetCurrentPSUnit(const G4String &unit)
virtual ~G4VScoringMesh()
virtual void DrawColumn(RunScore *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)=0
void RotateZ(G4double delta)
G4bool fGeometryHasBeenDestroyed
G4LogicalVolume * fMeshElementLogical
void SetCenterPosition(G4double centerPosition[3])
G4VPrimitiveScorer * fCurrentPS
void GetDivisionAxisNames(G4String divisionAxisNames[3])
G4ParallelWorldProcess * GetParallelWorldProcess() const
MeshScoreMap GetScoreMap() const
G4LogicalVolume * GetMeshElementLogical() const
MeshShape
G4String GetCurrentPSUnit()
G4String fDivisionAxisNames[3]
G4bool IsActive() const
double G4double
Definition: G4Types.hh:76
static constexpr double ps
Definition: G4SIunits.hh:172
void SetDrawPSName(const G4String &psname)
void SetNumberOfSegments(G4int nSegment[3])
std::map< G4String, RunScore * > MeshScoreMap
void SetVerboseLevel(G4int vl)
G4String fDrawPSName