source: trunk/source/digits_hits/utils/include/G4VScoringMesh.hh@ 1351

Last change on this file since 1351 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 8.4 KB
Line 
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,v 1.40 2010/07/27 01:44:54 akimura Exp $
28// GEANT4 tag $Name: $
29//
30
31#ifndef G4VScoringMesh_h
32#define G4VScoringMesh_h 1
33
34#include "globals.hh"
35#include "G4THitsMap.hh"
36#include "G4RotationMatrix.hh"
37
38class G4VPhysicalVolume;
39class G4MultiFunctionalDetector;
40class G4VPrimitiveScorer;
41class G4VSDFilter;
42class G4VScoreColorMap;
43
44#include <map>
45
46enum MeshShape { boxMesh, cylinderMesh, sphereMesh };
47typedef std::map<G4String,G4THitsMap<G4double>* > MeshScoreMap;
48// class description:
49//
50// This class represents a parallel world for interactive scoring purposes.
51//
52
53class G4VScoringMesh
54{
55 public:
56 G4VScoringMesh(G4String wName);
57 virtual ~G4VScoringMesh();
58
59 public: // with description
60 // a pure virtual function to construct this mesh geometry
61 virtual void Construct(G4VPhysicalVolume* fWorldPhys)=0;
62 // list infomration of this mesh
63 virtual void List() const;
64
65 public: // with description
66 // get the world name
67 inline const G4String& GetWorldName() const
68 { return fWorldName; }
69 // get whether this mesh is active or not
70 inline G4bool IsActive() const
71 { return fActive; }
72 // set an activity of this mesh
73 inline void Activate(G4bool vl = true)
74 { fActive = vl; }
75 // get the shape of this mesh
76 inline MeshShape GetShape() const
77 { return fShape; }
78 // accumulate hits in a registered primitive scorer
79 inline void Accumulate(G4THitsMap<G4double> * map);
80 // dump information of primitive socrers registered in this mesh
81 void Dump();
82 // draw a projected quantity on a current viewer
83 inline void DrawMesh(G4String psName,G4VScoreColorMap* colorMap,G4int axflg=111);
84 // draw a column of a quantity on a current viewer
85 inline void DrawMesh(G4String psName,G4int idxPlane,G4int iColumn,G4VScoreColorMap* colorMap);
86 // draw a projected quantity on a current viewer
87 virtual void Draw(std::map<G4int, G4double*> * map, G4VScoreColorMap* colorMap, G4int axflg=111) = 0;
88 // draw a column of a quantity on a current viewer
89 virtual void DrawColumn(std::map<G4int, G4double*> * map, G4VScoreColorMap* colorMap,
90 G4int idxProj, G4int idxColumn) = 0;
91 // reset registered primitive scorers
92 void ResetScore();
93
94 // set size of this mesh
95 void SetSize(G4double size[3]);
96 // get size of this mesh
97 G4ThreeVector GetSize() const;
98 // set position of center of this mesh
99 void SetCenterPosition(G4double centerPosition[3]);
100 // get position of center of this mesh
101 G4ThreeVector GetTranslation() const {return fCenterPosition;}
102 // set a rotation angle around the x axis
103 void RotateX(G4double delta);
104 // set a rotation angle around the y axis
105 void RotateY(G4double delta);
106 // set a rotation angle around the z axis
107 void RotateZ(G4double delta);
108 // get a rotation matrix
109 G4RotationMatrix GetRotationMatrix() const {
110 if(fRotationMatrix) return *fRotationMatrix;
111 else return G4RotationMatrix::IDENTITY;
112 }
113 // set number of segments of this mesh
114 void SetNumberOfSegments(G4int nSegment[3]);
115 // get number of segments of this mesh
116 void GetNumberOfSegments(G4int nSegment[3]);
117
118 // register a primitive scorer to the MFD & set it to the current primitive scorer
119 void SetPrimitiveScorer(G4VPrimitiveScorer * ps);
120 // register a filter to a current primtive scorer
121 void SetFilter(G4VSDFilter * filter);
122 // set a primitive scorer to the current one by the name
123 void SetCurrentPrimitiveScorer(G4String & name);
124 // find registered primitive scorer by the name
125 G4bool FindPrimitiveScorer(G4String & psname);
126 // get whether current primitive scorer is set or not
127 G4bool IsCurrentPrimitiveScorerNull() {
128 if(fCurrentPS == NULL) return true;
129 else return false;
130 }
131 // get unit of primitive scorer by the name
132 G4String GetPSUnit(G4String & psname);
133 // get unit of current primitive scorer
134 G4String GetCurrentPSUnit();
135 // set unit of current primitive scorer
136 void SetCurrentPSUnit(const G4String& unit);
137 // get unit value of primitive scorer by the name
138 G4double GetPSUnitValue(G4String & psname);
139 // set PS name to be drawn
140 void SetDrawPSName(G4String & psname) {fDrawPSName = psname;}
141
142 // get axis names of the hierarchical division in the divided order
143 void GetDivisionAxisNames(G4String divisionAxisNames[3]);
144
145 // set current primitive scorer to NULL
146 void SetNullToCurrentPrimitiveScorer() {fCurrentPS = NULL;}
147 // set verbose level
148 inline void SetVerboseLevel(G4int vl)
149 { verboseLevel = vl; }
150 // get the primitive scorer map
151 MeshScoreMap GetScoreMap() {return fMap;}
152 // get whether this mesh setup has been ready
153 inline G4bool ReadyForQuantity() const
154 { return (sizeIsSet && nMeshIsSet); }
155
156protected:
157 // get registered primitive socrer by the name
158 G4VPrimitiveScorer * GetPrimitiveScorer(G4String & name);
159
160protected:
161 G4String fWorldName;
162 G4VPrimitiveScorer * fCurrentPS;
163 G4bool fConstructed;
164 G4bool fActive;
165 MeshShape fShape;
166
167 G4double fSize[3];
168 G4ThreeVector fCenterPosition;
169 G4RotationMatrix * fRotationMatrix;
170 G4int fNSegment[3];
171
172 std::map<G4String, G4THitsMap<G4double>* > fMap;
173 G4MultiFunctionalDetector * fMFD;
174
175 G4int verboseLevel;
176
177 G4bool sizeIsSet;
178 G4bool nMeshIsSet;
179
180 G4String fDrawUnit;
181 G4double fDrawUnitValue;
182 G4String fDrawPSName;
183
184 G4String fDivisionAxisNames[3];
185};
186
187void G4VScoringMesh::Accumulate(G4THitsMap<G4double> * map)
188{
189 G4String psName = map->GetName();
190 std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
191 *(fMapItr->second) += *map;
192
193 if(verboseLevel > 9) {
194 G4cout << G4endl;
195 G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
196 G4cout << " PS name : " << psName << G4endl;
197 if(fMapItr == fMap.end()) {
198 G4cout << " "
199 << psName << " was not found." << G4endl;
200 } else {
201 G4cout << " map size : " << map->GetSize() << G4endl;
202 map->PrintAllHits();
203 }
204 G4cout << G4endl;
205 }
206}
207
208void G4VScoringMesh::DrawMesh(G4String psName,G4VScoreColorMap* colorMap,G4int axflg)
209{
210 fDrawPSName = psName;
211 std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
212 if(fMapItr!=fMap.end()) {
213 fDrawUnit = GetPSUnit(psName);
214 fDrawUnitValue = GetPSUnitValue(psName);
215 Draw(fMapItr->second->GetMap(), colorMap,axflg);
216 } else {
217 G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
218 }
219}
220
221void G4VScoringMesh::DrawMesh(G4String psName,G4int idxPlane,G4int iColumn,G4VScoreColorMap* colorMap)
222{
223 fDrawPSName = psName;
224 std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
225 if(fMapItr!=fMap.end()) {
226 fDrawUnit = GetPSUnit(psName);
227 fDrawUnitValue = GetPSUnitValue(psName);
228 DrawColumn(fMapItr->second->GetMap(),colorMap,idxPlane,iColumn);
229 } else {
230 G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
231 }
232}
233
234#endif
235
Note: See TracBrowser for help on using the repository browser.