Geant4 10.7.0
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 command-based scorer
50// For parallel world scorer, this class creates a parallel world mesh geometry
51//
52
54{
55public:
59 using MeshScoreMap = std::map< G4String, RunScore* >;
60
61public:
62 G4VScoringMesh(const G4String& wName);
63 virtual ~G4VScoringMesh();
64
65public: // with description
66 // a pure virtual function to construct this mesh geometry
67 void Construct(G4VPhysicalVolume* fWorldPhys);
68
69 void WorkerConstruct(G4VPhysicalVolume* fWorldPhys);
70
71protected:
72 virtual void SetupGeometry(G4VPhysicalVolume * fWorldPhys) = 0;
73
74public: // with description
75 // list infomration of this mesh
76 virtual void List() const;
77
78public: // with description
79 // get the world name
80 // If this ScoringMesh is for parallel world, it returns the name of the parallel world
81 // If this ScoringMesh is for real world logical volume, it returns name of logical volume
82 inline const G4String& GetWorldName() const
83 { return fWorldName; }
84 // get whether this mesh is active or not
85 inline G4bool IsActive() const
86 { return fActive; }
87 // set an activity of this mesh
88 inline void Activate(G4bool vl = true)
89 { fActive = vl; }
90 // get the shape of this mesh
91 inline MeshShape GetShape() const
92 { return fShape; }
93 // accumulate hits in a registered primitive scorer
96 // merge same kind of meshes
97 void Merge(const G4VScoringMesh * scMesh);
98 // dump information of primitive socrers registered in this mesh
99 void Dump();
100 // draw a projected quantity on a current viewer
101 void DrawMesh(const G4String& psName,G4VScoreColorMap* colorMap,G4int axflg=111);
102 // draw a column of a quantity on a current viewer
103 void DrawMesh(const G4String& psName,G4int idxPlane,G4int iColumn,G4VScoreColorMap* colorMap);
104 // draw a projected quantity on a current viewer
105 virtual void Draw(RunScore * map, G4VScoreColorMap* colorMap, G4int axflg=111) = 0;
106 // draw a column of a quantity on a current viewer
107 virtual void DrawColumn(RunScore * map, G4VScoreColorMap* colorMap,
108 G4int idxProj, G4int idxColumn) = 0;
109 // reset registered primitive scorers
110 void ResetScore();
111
112 // Following set/get methods make sense only for parallel world scoring mesh
113 // set size of this mesh
114 void SetSize(G4double size[3]);
115 // get size of this mesh
116 G4ThreeVector GetSize() const;
117 // set position of center of this mesh
118 void SetCenterPosition(G4double centerPosition[3]);
119 // get position of center of this mesh
121 // set a rotation angle around the x axis
122 void RotateX(G4double delta);
123 // set a rotation angle around the y axis
124 void RotateY(G4double delta);
125 // set a rotation angle around the z axis
126 void RotateZ(G4double delta);
127 // get a rotation matrix
130 else return G4RotationMatrix::IDENTITY;
131 }
132
133 // set number of segments of this mesh
134 void SetNumberOfSegments(G4int nSegment[3]);
135 // get number of segments of this mesh
136 void GetNumberOfSegments(G4int nSegment[3]);
137
138 // register a primitive scorer to the MFD & set it to the current primitive scorer
140 // register a filter to a current primtive scorer
141 void SetFilter(G4VSDFilter * filter);
142 // set a primitive scorer to the current one by the name
143 void SetCurrentPrimitiveScorer(const G4String & name);
144 // find registered primitive scorer by the name
145 G4bool FindPrimitiveScorer(const G4String & psname);
146 // get whether current primitive scorer is set or not
148 if(fCurrentPS == nullptr) return true;
149 else return false;
150 }
151 // get unit of primitive scorer by the name
152 G4String GetPSUnit(const G4String & psname);
153 // get unit of current primitive scorer
155 // set unit of current primitive scorer
156 void SetCurrentPSUnit(const G4String& unit);
157 // get unit value of primitive scorer by the name
158 G4double GetPSUnitValue(const G4String & psname);
159 // set PS name to be drawn
160 void SetDrawPSName(const G4String & psname) {fDrawPSName = psname;}
161
162 // get axis names of the hierarchical division in the divided order
163 void GetDivisionAxisNames(G4String divisionAxisNames[3]);
164
165 // set current primitive scorer to NULL
167 // set verbose level
168 inline void SetVerboseLevel(G4int vl)
169 { verboseLevel = vl; }
170 // get the primitive scorer map
171 inline MeshScoreMap GetScoreMap() const
172 { return fMap; }
173 // get whether this mesh setup has been ready
175 { return (sizeIsSet && nMeshIsSet); }
176
177//protected:
178 // get registered primitive socrer by the name
180
181protected:
187
192
195
197
200
204
206
208
209public:
211 { fMeshElementLogical = val; }
213 { return fMeshElementLogical; }
214
215protected:
218public:
220 { fParallelWorldProcess = proc; }
222 { return fParallelWorldProcess; }
224 {
226 fMeshElementLogical = nullptr;
227 }
228
229protected:
231public:
232 // Geometry hirarchy level (bottom = 0) to be used as the copy number
233 // This is used only for real-world scorer
234 inline void SetCopyNumberLevel(G4int val)
235 { copyNumberLevel = val; }
237 { return copyNumberLevel; }
238
239protected:
240 // This flag may be set to true for Probe scoring mesh.
241 // There is no public set method for this boolean flag, but it should be set to true
242 // through SetMaterial() method of Probe scoring mesh.
244public:
246 { return layeredMassFlg; }
247};
248
249#endif
250
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
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
void WorkerConstruct(G4VPhysicalVolume *fWorldPhys)
G4String fDrawPSName
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 SetCurrentPSUnit(const G4String &unit)
void SetParallelWorldProcess(G4ParallelWorldProcess *proc)
G4MultiFunctionalDetector * fMFD
void SetCurrentPrimitiveScorer(const G4String &name)
const G4String & GetWorldName() const
std::map< G4String, RunScore * > MeshScoreMap
G4String fDivisionAxisNames[3]
G4LogicalVolume * fMeshElementLogical
void SetPrimitiveScorer(G4VPrimitiveScorer *ps)
G4VPrimitiveScorer * fCurrentPS
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]
G4ParallelWorldProcess * fParallelWorldProcess
MeshScoreMap GetScoreMap() const
void Merge(const G4VScoringMesh *scMesh)
virtual ~G4VScoringMesh()
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)
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