Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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$
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
40class G4VSDFilter;
42
43#include <map>
44
46typedef std::map<G4String,G4THitsMap<G4double>* > MeshScoreMap;
47// class description:
48//
49// This class represents a parallel world for interactive scoring purposes.
50//
51
53{
54 public:
55 G4VScoringMesh(const G4String& wName);
56 virtual ~G4VScoringMesh();
57
58 public: // with description
59 // a pure virtual function to construct this mesh geometry
60 virtual void Construct(G4VPhysicalVolume* fWorldPhys)=0;
61 // list infomration of this mesh
62 virtual void List() const;
63
64 public: // with description
65 // get the world name
66 inline const G4String& GetWorldName() const
67 { return fWorldName; }
68 // get whether this mesh is active or not
69 inline G4bool IsActive() const
70 { return fActive; }
71 // set an activity of this mesh
72 inline void Activate(G4bool vl = true)
73 { fActive = vl; }
74 // get the shape of this mesh
75 inline MeshShape GetShape() const
76 { return fShape; }
77 // accumulate hits in a registered primitive scorer
78 inline void Accumulate(G4THitsMap<G4double> * map);
79 // dump information of primitive socrers registered in this mesh
80 void Dump();
81 // draw a projected quantity on a current viewer
82 void DrawMesh(const G4String& psName,G4VScoreColorMap* colorMap,G4int axflg=111);
83 // draw a column of a quantity on a current viewer
84 void DrawMesh(const G4String& psName,G4int idxPlane,G4int iColumn,G4VScoreColorMap* colorMap);
85 // draw a projected quantity on a current viewer
86 virtual void Draw(std::map<G4int, G4double*> * map, G4VScoreColorMap* colorMap, G4int axflg=111) = 0;
87 // draw a column of a quantity on a current viewer
88 virtual void DrawColumn(std::map<G4int, G4double*> * map, G4VScoreColorMap* colorMap,
89 G4int idxProj, G4int idxColumn) = 0;
90 // reset registered primitive scorers
91 void ResetScore();
92
93 // set size of this mesh
94 void SetSize(G4double size[3]);
95 // get size of this mesh
96 G4ThreeVector GetSize() const;
97 // set position of center of this mesh
98 void SetCenterPosition(G4double centerPosition[3]);
99 // get position of center of this mesh
101 // set a rotation angle around the x axis
102 void RotateX(G4double delta);
103 // set a rotation angle around the y axis
104 void RotateY(G4double delta);
105 // set a rotation angle around the z axis
106 void RotateZ(G4double delta);
107 // get a rotation matrix
110 else return G4RotationMatrix::IDENTITY;
111 }
112 // set number of segments of this mesh
113 void SetNumberOfSegments(G4int nSegment[3]);
114 // get number of segments of this mesh
115 void GetNumberOfSegments(G4int nSegment[3]);
116
117 // register a primitive scorer to the MFD & set it to the current primitive scorer
119 // register a filter to a current primtive scorer
120 void SetFilter(G4VSDFilter * filter);
121 // set a primitive scorer to the current one by the name
122 void SetCurrentPrimitiveScorer(const G4String & name);
123 // find registered primitive scorer by the name
124 G4bool FindPrimitiveScorer(const G4String & psname);
125 // get whether current primitive scorer is set or not
127 if(fCurrentPS == NULL) return true;
128 else return false;
129 }
130 // get unit of primitive scorer by the name
131 G4String GetPSUnit(const G4String & psname);
132 // get unit of current primitive scorer
134 // set unit of current primitive scorer
135 void SetCurrentPSUnit(const G4String& unit);
136 // get unit value of primitive scorer by the name
137 G4double GetPSUnitValue(const G4String & psname);
138 // set PS name to be drawn
139 void SetDrawPSName(const G4String & psname) {fDrawPSName = psname;}
140
141 // get axis names of the hierarchical division in the divided order
142 void GetDivisionAxisNames(G4String divisionAxisNames[3]);
143
144 // set current primitive scorer to NULL
146 // set verbose level
147 inline void SetVerboseLevel(G4int vl)
148 { verboseLevel = vl; }
149 // get the primitive scorer map
151 // get whether this mesh setup has been ready
153 { return (sizeIsSet && nMeshIsSet); }
154
155protected:
156 // get registered primitive socrer by the name
158
159protected:
165
170
171 std::map<G4String, G4THitsMap<G4double>* > fMap;
173
175
178
182
184};
185
187{
188 G4String psName = map->GetName();
189 std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
190 *(fMapItr->second) += *map;
191
192 if(verboseLevel > 9) {
193 G4cout << G4endl;
194 G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
195 G4cout << " PS name : " << psName << G4endl;
196 if(fMapItr == fMap.end()) {
197 G4cout << " "
198 << psName << " was not found." << G4endl;
199 } else {
200 G4cout << " map size : " << map->GetSize() << G4endl;
201 map->PrintAllHits();
202 }
203 G4cout << G4endl;
204 }
205}
206
207#endif
208
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
std::map< G4String, G4THitsMap< G4double > * > MeshScoreMap
MeshShape
@ sphereMesh
@ boxMesh
@ cylinderMesh
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static DLL_API const HepRotation IDENTITY
Definition: Rotation.h:369
virtual size_t GetSize() const
Definition: G4THitsMap.hh:86
virtual void PrintAllHits()
Definition: G4THitsMap.hh:193
void SetVerboseLevel(G4int vl)
void Activate(G4bool vl=true)
void SetFilter(G4VSDFilter *filter)
std::map< G4String, G4THitsMap< G4double > * > fMap
G4bool ReadyForQuantity() const
void SetNullToCurrentPrimitiveScorer()
G4RotationMatrix * fRotationMatrix
virtual void Construct(G4VPhysicalVolume *fWorldPhys)=0
G4ThreeVector GetTranslation() const
G4ThreeVector GetSize() const
MeshShape GetShape() const
virtual void List() const
G4double fDrawUnitValue
G4String fDrawPSName
G4bool IsActive() const
void SetDrawPSName(const G4String &psname)
void RotateY(G4double delta)
void GetNumberOfSegments(G4int nSegment[3])
G4RotationMatrix GetRotationMatrix() const
G4double GetPSUnitValue(const G4String &psname)
G4String GetPSUnit(const G4String &psname)
virtual void Draw(std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0
void SetCurrentPSUnit(const G4String &unit)
MeshScoreMap GetScoreMap()
G4MultiFunctionalDetector * fMFD
void SetCurrentPrimitiveScorer(const G4String &name)
const G4String & GetWorldName() const
G4String fDivisionAxisNames[3]
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]
virtual ~G4VScoringMesh()
virtual void DrawColumn(std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)=0
G4bool IsCurrentPrimitiveScorerNull()
void SetCenterPosition(G4double centerPosition[3])
void RotateX(G4double delta)
void SetSize(G4double size[3])
G4bool FindPrimitiveScorer(const G4String &psname)
void RotateZ(G4double delta)
G4VPrimitiveScorer * GetPrimitiveScorer(const G4String &name)
G4ThreeVector fCenterPosition