Geant4  10.01.p02
G4VScoringMesh.cc
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: G4VScoringMesh.cc 83518 2014-08-27 12:57:10Z gcosmo $
28 //
29 // ---------------------------------------------------------------------
30 // Modifications
31 // 17-Apr-2012 T.Aso SetSize() and SetNumberOfSegments() is not allowed
32 // to call twice in same geometrical mesh. Add warning
33 // message to notify.
34 //
35 // ---------------------------------------------------------------------
36 
37 #include "G4VScoringMesh.hh"
38 
39 #include "G4SystemOfUnits.hh"
40 #include "G4VPhysicalVolume.hh"
42 #include "G4VPrimitiveScorer.hh"
43 #include "G4VSDFilter.hh"
44 #include "G4SDManager.hh"
45 
47  : fWorldName(wName),fCurrentPS(0),fConstructed(false),fActive(true),fShape(undefinedMesh),
48  fRotationMatrix(0), fMFD(new G4MultiFunctionalDetector(wName)),
49  verboseLevel(0),sizeIsSet(false),nMeshIsSet(false),
50  fDrawUnit(""), fDrawUnitValue(1.), fMeshElementLogical(0),
51  fParallelWorldProcess(0), fGeometryHasBeenDestroyed(false)
52 {
54 
55  fSize[0] = fSize[1] = fSize[2] = 0.;
56  fNSegment[0] = fNSegment[1] = fNSegment[2] = 1;
58 }
59 
61 {
62  ;
63 }
64 
66  if(verboseLevel > 9) G4cout << "G4VScoringMesh::ResetScore() is called." << G4endl;
67  std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.begin();
68  for(; itr != fMap.end(); itr++) {
69  if(verboseLevel > 9) G4cout << "G4VScoringMesh::ResetScore()" << itr->first << G4endl;
70  itr->second->clear();
71  }
72 }
74  if ( !sizeIsSet ){
75  for(int i = 0; i < 3; i++) fSize[i] = size[i];
76  sizeIsSet = true;
77  }else{
78  G4String message = " The size of scoring mesh can not be changed.";
79  G4Exception("G4VScoringMesh::SetSize()",
80  "DigiHitsUtilsScoreVScoringMesh000", JustWarning,
81  message);
82  }
83 }
85  if(sizeIsSet)
86  return G4ThreeVector(fSize[0], fSize[1], fSize[2]);
87  else
88  return G4ThreeVector(0., 0., 0.);
89 }
91  fCenterPosition = G4ThreeVector(centerPosition[0], centerPosition[1], centerPosition[2]);
92 }
94  if ( !nMeshIsSet ){
95  for(int i = 0; i < 3; i++) fNSegment[i] = nSegment[i];
96  nMeshIsSet = true;
97  } else {
98  G4String message = " The size of scoring segments can not be changed.";
99  G4Exception("G4VScoringMesh::SetNumberOfSegments()",
100  "DigiHitsUtilsScoreVScoringMesh000", JustWarning,
101  message);
102  }
103 }
105  for(int i = 0; i < 3; i++) nSegment[i] = fNSegment[i];
106 }
109  fRotationMatrix->rotateX(delta);
110 }
111 
114  fRotationMatrix->rotateY(delta);
115 }
116 
119  fRotationMatrix->rotateZ(delta);
120 }
121 
123 
124  if(!ReadyForQuantity())
125  {
126  G4cerr << "ERROR : G4VScoringMesh::SetPrimitiveScorer() : " << ps->GetName()
127  << " does not yet have mesh size or number of bins. Set them first." << G4endl
128  << "This Method is ignored." << G4endl;
129  return;
130  }
131  if(verboseLevel > 0) G4cout << "G4VScoringMesh::SetPrimitiveScorer() : "
132  << ps->GetName() << " is registered."
133  << " 3D size: ("
134  << fNSegment[0] << ", "
135  << fNSegment[1] << ", "
136  << fNSegment[2] << ")" << G4endl;
137 
138  ps->SetNijk(fNSegment[0], fNSegment[1], fNSegment[2]);
139  fCurrentPS = ps;
140  fMFD->RegisterPrimitive(ps);
142  fMap[ps->GetName()] = map;
143 }
144 
146 
147  if(fCurrentPS == 0) {
148  G4cerr << "ERROR : G4VScoringMesh::SetSDFilter() : a quantity must be defined first. This method is ignored." << G4endl;
149  return;
150  }
151  if(verboseLevel > 0) G4cout << "G4VScoringMesh::SetFilter() : "
152  << filter->GetName()
153  << " is set to "
154  << fCurrentPS->GetName() << G4endl;
155 
156  G4VSDFilter* oldFilter = fCurrentPS->GetFilter();
157  if(oldFilter)
158  {
159  G4cout << "WARNING : G4VScoringMesh::SetFilter() : " << oldFilter->GetName()
160  << " is overwritten by " << filter->GetName() << G4endl;
161  }
162  fCurrentPS->SetFilter(filter);
163 }
164 
167  if(fCurrentPS == 0) {
168  G4cerr << "ERROR : G4VScoringMesh::SetCurrentPrimitiveScorer() : The primitive scorer <"
169  << name << "> does not found." << G4endl;
170  }
171 }
172 
174  std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.find(psname);;
175  if(itr == fMap.end()) return false;
176  return true;
177 }
178 
180  std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.find(psname);;
181  if(itr == fMap.end()) {
182  return G4String("");
183  } else {
184  return GetPrimitiveScorer(psname)->GetUnit();
185  }
186 }
187 
189  G4String unit = "";
190  if(fCurrentPS == 0) {
191  G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
192  msg += " Current primitive scorer is null.";
193  G4cerr << msg << G4endl;
194  }else{
195  unit = fCurrentPS->GetUnit();
196  }
197  return unit;
198 }
199 
201  if(fCurrentPS == 0) {
202  G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
203  msg += " Current primitive scorer is null.";
204  G4cerr << msg << G4endl;
205  }else{
206  fCurrentPS->SetUnit(unit);
207  }
208 }
209 
211  std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.find(psname);;
212  if(itr == fMap.end()) {
213  return 1.;
214  } else {
215  return GetPrimitiveScorer(psname)->GetUnitValue();
216  }
217 }
218 
219 void G4VScoringMesh::GetDivisionAxisNames(G4String divisionAxisNames[3]) {
220  for(int i = 0; i < 3; i++) divisionAxisNames[i] = fDivisionAxisNames[i];
221 }
222 
224  if(fMFD == 0) return 0;
225 
227  for(G4int i = 0; i < nps; i++) {
229  if(name == ps->GetName()) return ps;
230  }
231 
232  return 0;
233 }
234 void G4VScoringMesh::List() const {
235 
236  G4cout << " # of segments: ("
237  << fNSegment[0] << ", "
238  << fNSegment[1] << ", "
239  << fNSegment[2] << ")"
240  << G4endl;
241  G4cout << " displacement: ("
242  << fCenterPosition.x()/cm << ", "
243  << fCenterPosition.y()/cm << ", "
244  << fCenterPosition.z()/cm << ") [cm]"
245  << G4endl;
246  if(fRotationMatrix != 0) {
247  G4cout << " rotation matrix: "
248  << fRotationMatrix->xx() << " "
249  << fRotationMatrix->xy() << " "
250  << fRotationMatrix->xz() << G4endl
251  << " "
252  << fRotationMatrix->yx() << " "
253  << fRotationMatrix->yy() << " "
254  << fRotationMatrix->yz() << G4endl
255  << " "
256  << fRotationMatrix->zx() << " "
257  << fRotationMatrix->zy() << " "
258  << fRotationMatrix->zz() << G4endl;
259  }
260 
261 
262  G4cout << " registered primitve scorers : " << G4endl;
264  G4VPrimitiveScorer * ps;
265  for(int i = 0; i < nps; i++) {
266  ps = fMFD->GetPrimitive(i);
267  G4cout << " " << i << " " << ps->GetName();
268  if(ps->GetFilter() != 0) G4cout << " with " << ps->GetFilter()->GetName();
269  G4cout << G4endl;
270  }
271 
272 
273 }
274 
276  G4cout << "scoring mesh name: " << fWorldName << G4endl;
277 
278  G4cout << "# of G4THitsMap : " << fMap.size() << G4endl;
279  std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.begin();
280  for(; itr != fMap.end(); itr++) {
281  G4cout << "[" << itr->first << "]" << G4endl;
282  itr->second->PrintAllHits();
283  }
284  G4cout << G4endl;
285 
286 }
287 
288 
289 void G4VScoringMesh::DrawMesh(const G4String& psName,G4VScoreColorMap* colorMap,G4int axflg)
290 {
291  fDrawPSName = psName;
292  std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
293  if(fMapItr!=fMap.end()) {
294  fDrawUnit = GetPSUnit(psName);
295  fDrawUnitValue = GetPSUnitValue(psName);
296  Draw(fMapItr->second->GetMap(), colorMap,axflg);
297  } else {
298  G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
299  }
300 }
301 
302 void G4VScoringMesh::DrawMesh(const G4String& psName,G4int idxPlane,G4int iColumn,G4VScoreColorMap* colorMap)
303 {
304  fDrawPSName = psName;
305  std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
306  if(fMapItr!=fMap.end()) {
307  fDrawUnit = GetPSUnit(psName);
308  fDrawUnitValue = GetPSUnitValue(psName);
309  DrawColumn(fMapItr->second->GetMap(),colorMap,idxPlane,iColumn);
310  } else {
311  G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
312  }
313 }
314 
316 {
317  G4String psName = map->GetName();
318  std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
319  *(fMapItr->second) += *map;
320 
321  if(verboseLevel > 9) {
322  G4cout << G4endl;
323  G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
324  G4cout << " PS name : " << psName << G4endl;
325  if(fMapItr == fMap.end()) {
326  G4cout << " "
327  << psName << " was not found." << G4endl;
328  } else {
329  G4cout << " map size : " << map->GetSize() << G4endl;
330  map->PrintAllHits();
331  }
332  G4cout << G4endl;
333  }
334 }
335 
337 {
338  if(fConstructed) {
340  SetupGeometry(fWorldPhys);
342  }
343  if(verboseLevel > 0)
344  G4cout << fWorldPhys->GetName() << " --- All quantities are reset."
345  << G4endl;
346  ResetScore();
347  }
348  else {
349  fConstructed = true;
350  SetupGeometry(fWorldPhys);
351  }
352 }
353 
355 {
356  if(fConstructed) {
360  }
361 
362  if(verboseLevel > 0)
363  G4cout << fWorldPhys->GetName() << " --- All quantities are reset." << G4endl;
364  ResetScore();
365 
366  } else {
367  fConstructed = true;
369  }
370 }
371 
373 {
374  const MeshScoreMap scMap = scMesh->GetScoreMap();
375 
376  MeshScoreMap::const_iterator fMapItr = fMap.begin();
377  MeshScoreMap::const_iterator mapItr = scMap.begin();
378  for(; fMapItr != fMap.end(); fMapItr++) {
379  if(verboseLevel > 9) G4cout << "G4VScoringMesh::Merge()" << fMapItr->first << G4endl;
380  *(fMapItr->second) += *(mapItr->second);
381  mapItr++;
382  }
383 }
384 
static const double cm
Definition: G4SIunits.hh:106
void GetNumberOfSegments(G4int nSegment[3])
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
G4VScoringMesh(const G4String &wName)
void DrawMesh(const G4String &psName, G4VScoreColorMap *colorMap, G4int axflg=111)
G4ThreeVector fCenterPosition
G4String GetName() const
void Merge(const G4VScoringMesh *scMesh)
G4bool FindPrimitiveScorer(const G4String &psname)
std::map< G4String, G4THitsMap< G4double > * > MeshScoreMap
CLHEP::Hep3Vector G4ThreeVector
void RotateY(G4double delta)
G4VPrimitiveScorer * GetPrimitiveScorer(const G4String &name)
CLHEP::HepRotation G4RotationMatrix
void SetPrimitiveScorer(G4VPrimitiveScorer *ps)
void SetCurrentPrimitiveScorer(const G4String &name)
G4String name
Definition: TRTMaterials.hh:40
G4MultiFunctionalDetector * fMFD
void SetSize(G4double size[3])
G4double GetPSUnitValue(const G4String &psname)
void SetFilter(G4VSDFilter *f)
void Accumulate(G4THitsMap< G4double > *map)
void RotateX(G4double delta)
G4double fSize[3]
G4String GetName() const
Definition: G4VSDFilter.hh:57
void SetFilter(G4VSDFilter *filter)
void WorkerConstruct(G4VPhysicalVolume *fWorldPhys)
G4ThreeVector GetSize() const
G4String GetPSUnit(const G4String &psname)
int G4int
Definition: G4Types.hh:78
void Construct(G4VPhysicalVolume *fWorldPhys)
void SetNijk(G4int i, G4int j, G4int k)
G4double fDrawUnitValue
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
virtual void List() const
G4VPrimitiveScorer * GetPrimitive(G4int id) const
G4double GetUnitValue() const
G4RotationMatrix * fRotationMatrix
bool G4bool
Definition: G4Types.hh:79
virtual size_t GetSize() const
Definition: G4THitsMap.hh:86
virtual void DrawColumn(std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)=0
G4bool ReadyForQuantity() const
virtual void SetupGeometry(G4VPhysicalVolume *fWorldPhys)=0
void SetCurrentPSUnit(const G4String &unit)
virtual ~G4VScoringMesh()
void RotateZ(G4double delta)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4bool fGeometryHasBeenDestroyed
G4LogicalVolume * fMeshElementLogical
void AddNewDetector(G4VSensitiveDetector *aSD)
Definition: G4SDManager.cc:67
virtual void Draw(std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0
void SetCenterPosition(G4double centerPosition[3])
G4VPrimitiveScorer * fCurrentPS
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
void GetDivisionAxisNames(G4String divisionAxisNames[3])
virtual void PrintAllHits()
Definition: G4THitsMap.hh:193
MeshScoreMap GetScoreMap() const
#define G4endl
Definition: G4ios.hh:61
G4String GetCurrentPSUnit()
G4String fDivisionAxisNames[3]
double G4double
Definition: G4Types.hh:76
const G4String & GetUnit() const
void SetUnit(const G4String &unit)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
void SetNumberOfSegments(G4int nSegment[3])
std::map< G4String, G4THitsMap< G4double > * > fMap
G4String fDrawPSName
G4GLOB_DLL std::ostream G4cerr
G4VSDFilter * GetFilter() const