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

Last change on this file since 1331 was 1228, checked in by garnier, 14 years ago

update geant4.9.3 tag

File size: 7.7 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.29 2009/10/12 04:11:25 akimura Exp $
28// GEANT4 tag $Name: geant4-09-03 $
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  ~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  // set positions to segment this mesh
118  inline void SetSegmentPositions(std::vector<G4double> & sp) {fSegmentPositions = sp;}
119
120  // register a primitive scorer to the MFD & set it to the current primitive scorer
121  void SetPrimitiveScorer(G4VPrimitiveScorer * ps);
122  // register a filter to a current primtive scorer
123  void SetFilter(G4VSDFilter * filter);
124  // set a primitive scorer to the current one by the name
125  void SetCurrentPrimitiveScorer(G4String & name);
126  // find registered primitive scorer by the name
127  G4bool FindPrimitiveScorer(G4String & psname);
128  // get whether current primitive scorer is set or not
129  G4bool IsCurrentPrimitiveScorerNull() {
130    if(fCurrentPS == NULL) return true;
131    else return false;
132  }
133  // set current  primitive scorer to NULL
134  void SetNullToCurrentPrimitiveScorer() {fCurrentPS = NULL;}
135  // set verbose level
136  inline void SetVerboseLevel(G4int vl) 
137  { verboseLevel = vl; }
138  // get the primitive scorer map
139  MeshScoreMap GetScoreMap() {return fMap;}
140  // get whether this mesh setup has been ready
141  inline G4bool ReadyForQuantity() const
142  { return (sizeIsSet && nMeshIsSet); }
143
144protected:
145  // get registered primitive socrer by the name
146  G4VPrimitiveScorer * GetPrimitiveScorer(G4String & name);
147
148protected:
149  G4String  fWorldName;
150  G4VPrimitiveScorer * fCurrentPS;
151  G4bool    fConstructed;
152  G4bool    fActive;
153  MeshShape fShape;
154
155  G4double fSize[3];
156  G4ThreeVector fCenterPosition;
157  G4RotationMatrix * fRotationMatrix;
158  G4int fNSegment[3];
159  std::vector<G4double> fSegmentPositions;
160
161  std::map<G4String, G4THitsMap<G4double>* > fMap;
162  G4MultiFunctionalDetector * fMFD;
163
164  G4int verboseLevel;
165
166  G4bool sizeIsSet;
167  G4bool nMeshIsSet;
168
169};
170
171void G4VScoringMesh::Accumulate(G4THitsMap<G4double> * map)
172{
173  G4String psName = map->GetName();
174  std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
175  *(fMapItr->second) += *map;
176
177  if(verboseLevel > 9) {
178    G4cout << G4endl;
179    G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
180    G4cout << "  PS name : " << psName << G4endl;
181    if(fMapItr == fMap.end()) {
182      G4cout << "  "
183             << psName << " was not found." << G4endl;
184    } else {
185      G4cout << "  map size : " << map->GetSize() << G4endl;
186      map->PrintAllHits();
187    }
188    G4cout << G4endl;
189  }
190}
191
192void G4VScoringMesh::DrawMesh(G4String psName,G4VScoreColorMap* colorMap,G4int axflg)
193{
194  std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
195  if(fMapItr!=fMap.end())
196  { Draw(fMapItr->second->GetMap(),colorMap,axflg); }
197  else
198  { G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl; }
199}
200
201void G4VScoringMesh::DrawMesh(G4String psName,G4int idxPlane,G4int iColumn,G4VScoreColorMap* colorMap)
202{
203  std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
204  if(fMapItr!=fMap.end())
205  { DrawColumn(fMapItr->second->GetMap(),colorMap,idxPlane,iColumn); }
206  else
207  { G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl; }
208}
209
210#endif
211
Note: See TracBrowser for help on using the repository browser.