Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
28
29#ifndef G4VScoringMesh_h
30#define G4VScoringMesh_h 1
31
32#include "globals.hh"
33#include "G4THitsMap.hh"
34#include "G4RotationMatrix.hh"
35#include "G4StatDouble.hh"
36
38class G4LogicalVolume;
41class G4VSDFilter;
44
45#include <map>
46
47// class description:
48//
49// This class represents a multi-functional detector to be used by
50// command-based scorer For parallel world scorer, this class creates a
51// parallel world mesh geometry
52//
53
55{
56 public:
57 enum class MeshShape
58 {
59 box,
61 sphere,
63 probe,
64 undefined = -1
65 };
68 using MeshScoreMap = std::map<G4String, RunScore*>;
69
70 public:
71 G4VScoringMesh(const G4String& wName);
72 virtual ~G4VScoringMesh() = default;
73
74 public: // with description
75 virtual void Construct(G4VPhysicalVolume* fWorldPhys);
76 virtual void WorkerConstruct(G4VPhysicalVolume* fWorldPhys);
77
78 // list infomration of this mesh
79 virtual void List() const;
80
81 // get the world name
82 // If this ScoringMesh is for parallel world, it returns the name of the
83 // parallel world If this ScoringMesh is for real world logical volume, it
84 // returns name of logical volume
85 inline const G4String& GetWorldName() const { return fWorldName; }
86 // get whether this mesh is active or not
87 inline G4bool IsActive() const { return fActive; }
88 // set an activity of this mesh
89 inline void Activate(G4bool vl = true) { fActive = vl; }
90 // get the shape of this mesh
91 inline MeshShape GetShape() const { return fShape; }
92 // accumulate hits in a registered primitive scorer
95 // merge same kind of meshes
96 void Merge(const G4VScoringMesh* scMesh);
97 // dump information of primitive socrers registered in this mesh
98 void Dump();
99 // draw a projected quantity on a current viewer
100 void DrawMesh(const G4String& psName, G4VScoreColorMap* colorMap,
101 G4int axflg = 111);
102 // draw a column of a quantity on a current viewer
103 void DrawMesh(const G4String& psName, G4int idxPlane, G4int iColumn,
104 G4VScoreColorMap* colorMap);
105 // draw a projected quantity on a current viewer
106 virtual void Draw(RunScore* map, G4VScoreColorMap* colorMap,
107 G4int axflg = 111) = 0;
108 // draw a column of a quantity on a current viewer
109 virtual void DrawColumn(RunScore* map, G4VScoreColorMap* colorMap,
110 G4int idxProj, G4int idxColumn) = 0;
111 // reset registered primitive scorers
112 void ResetScore();
113
114 // Following set/get methods make sense only for parallel world scoring mesh
115 // set size of this mesh
116 void SetSize(G4double size[3]);
117 // get size of this mesh
118 G4ThreeVector GetSize() const;
119 // set starting and span angles (used only for tube segment)
121 // get angles (used only for tube segment)
122 inline G4double GetStartAngle() const { return fAngle[0]; }
123 inline G4double GetAngleSpan() const { return fAngle[1]; }
124 // set position of center of this mesh
125 void SetCenterPosition(G4double centerPosition[3]);
126 // get position of center of this mesh
128 // set a rotation angle around the x axis
129 void RotateX(G4double delta);
130 // set a rotation angle around the y axis
131 void RotateY(G4double delta);
132 // set a rotation angle around the z axis
133 void RotateZ(G4double delta);
134 // get a rotation matrix
136 {
137 if(fRotationMatrix != nullptr)
138 return *fRotationMatrix;
140 }
141
142 // set number of segments of this mesh
143 void SetNumberOfSegments(G4int nSegment[3]);
144 // get number of segments of this mesh
145 void GetNumberOfSegments(G4int nSegment[3]);
146
147 // register a primitive scorer to the MFD & set it to the current primitive
148 // scorer
150 // register a filter to a current primtive scorer
151 void SetFilter(G4VSDFilter* filter);
152 // set a primitive scorer to the current one by the name
153 void SetCurrentPrimitiveScorer(const G4String& name);
154 // find registered primitive scorer by the name
155 G4bool FindPrimitiveScorer(const G4String& psname);
156 // get whether current primitive scorer is set or not
158 {
159 return fCurrentPS == nullptr;
160 }
161 // get unit of primitive scorer by the name
162 G4String GetPSUnit(const G4String& psname);
163 // get unit of current primitive scorer
165 // set unit of current primitive scorer
166 void SetCurrentPSUnit(const G4String& unit);
167 // get unit value of primitive scorer by the name
168 G4double GetPSUnitValue(const G4String& psname);
169 // set PS name to be drawn
170 void SetDrawPSName(const G4String& psname) { fDrawPSName = psname; }
171
172 // get axis names of the hierarchical division in the divided order
173 void GetDivisionAxisNames(G4String divisionAxisNames[3]);
174
175 // set current primitive scorer to NULL
177 // set verbose level
178 inline void SetVerboseLevel(G4int vl) { verboseLevel = vl; }
179 // get the primitive scorer map
180 inline MeshScoreMap GetScoreMap() const { return fMap; }
181 // get whether this mesh setup has been ready
182 inline G4bool ReadyForQuantity() const { return (sizeIsSet && nMeshIsSet); }
183
184 // protected:
185 // get registered primitive socrer by the name
187
189 {
191 }
193 {
194 return fMeshElementLogical;
195 }
196
198 {
200 }
206 {
208 fMeshElementLogical = nullptr;
209 }
210
211 // Geometry hirarchy level (bottom = 0) to be used as the copy number
212 // This is used only for real-world scorer
213 inline void SetCopyNumberLevel(G4int val) { copyNumberLevel = val; }
214 inline G4int GetCopyNumberLevel() const { return copyNumberLevel; }
215
217
218 protected:
219 // a pure virtual function to construct this mesh geometry
220 virtual void SetupGeometry(G4VPhysicalVolume* fWorldPhys) = 0;
221
222 protected:
228
234
237
239
242
246
248
250
253
255
256 // This flag may be set to true for Probe scoring mesh.
257 // There is no public set method for this boolean flag, but it should be set
258 // to true through SetMaterial() method of Probe scoring mesh.
260};
261
262#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
static DLL_API const HepRotation IDENTITY
Definition Rotation.h:366
virtual ~G4VScoringMesh()=default
void SetVerboseLevel(G4int vl)
void Activate(G4bool vl=true)
void SetFilter(G4VSDFilter *filter)
G4bool ReadyForQuantity() const
void SetNullToCurrentPrimitiveScorer()
G4RotationMatrix * fRotationMatrix
G4ThreeVector GetTranslation() const
G4ThreeVector GetSize() const
MeshShape GetShape() const
virtual void SetupGeometry(G4VPhysicalVolume *fWorldPhys)=0
virtual void List() const
G4LogicalVolume * GetMeshElementLogical() const
G4double fDrawUnitValue
virtual void WorkerConstruct(G4VPhysicalVolume *fWorldPhys)
G4double GetStartAngle() const
G4bool IsActive() const
void SetDrawPSName(const G4String &psname)
G4bool fGeometryHasBeenDestroyed
MeshScoreMap fMap
void GeometryHasBeenDestroyed()
void RotateY(G4double delta)
void GetNumberOfSegments(G4int nSegment[3])
G4RotationMatrix GetRotationMatrix() const
G4int GetCopyNumberLevel() const
G4double GetPSUnitValue(const G4String &psname)
G4ParallelWorldProcess * GetParallelWorldProcess() const
G4String GetPSUnit(const G4String &psname)
void SetAngles(G4double, G4double)
void SetCurrentPSUnit(const G4String &unit)
void SetParallelWorldProcess(G4ParallelWorldProcess *proc)
G4MultiFunctionalDetector * fMFD
void SetCurrentPrimitiveScorer(const G4String &name)
const G4String & GetWorldName() const
G4String fDivisionAxisNames[3]
G4LogicalVolume * fMeshElementLogical
void SetPrimitiveScorer(G4VPrimitiveScorer *ps)
std::map< G4String, RunScore * > MeshScoreMap
G4VPrimitiveScorer * fCurrentPS
G4double fAngle[2]
G4String GetCurrentPSUnit()
void GetDivisionAxisNames(G4String divisionAxisNames[3])
void DrawMesh(const G4String &psName, G4VScoreColorMap *colorMap, G4int axflg=111)
void Accumulate(G4THitsMap< G4double > *map)
void SetNumberOfSegments(G4int nSegment[3])
G4double fSize[3]
G4VScoringMesh(const G4String &wName)
G4ParallelWorldProcess * fParallelWorldProcess
G4double GetAngleSpan() const
MeshScoreMap GetScoreMap() const
void Merge(const G4VScoringMesh *scMesh)
virtual void Draw(RunScore *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0
G4bool IsCurrentPrimitiveScorerNull()
void SetCenterPosition(G4double centerPosition[3])
void SetCopyNumberLevel(G4int val)
void RotateX(G4double delta)
virtual void Construct(G4VPhysicalVolume *fWorldPhys)
void SetSize(G4double size[3])
virtual void DrawColumn(RunScore *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)=0
G4bool FindPrimitiveScorer(const G4String &psname)
void SetMeshElementLogical(G4LogicalVolume *val)
void RotateZ(G4double delta)
G4bool LayeredMassFlg()
G4VPrimitiveScorer * GetPrimitiveScorer(const G4String &name)
G4ThreeVector fCenterPosition