Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4VScoringMesh.cc
이 파일의 문서화 페이지로 가기
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 99262 2016-09-09 13:18:19Z 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 #include "G4THitsMap.hh"
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(nullptr),fConstructed(false),fActive(true),fShape(undefinedMesh),
48  fRotationMatrix(nullptr), fMFD(new G4MultiFunctionalDetector(wName)),
49  verboseLevel(0),sizeIsSet(false),nMeshIsSet(false),
50  fDrawUnit(""), fDrawUnitValue(1.), fMeshElementLogical(nullptr),
51  fParallelWorldProcess(nullptr), 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 
65  if(verboseLevel > 9) G4cout << "G4VScoringMesh::ResetScore() is called." << G4endl;
66  for(auto mp : fMap)
67  {
68  if(verboseLevel > 9) G4cout << "G4VScoringMesh::ResetScore()" << mp.first << G4endl;
69  mp.second->clear();
70  }
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() : "
127  << prs->GetName()
128  << " does not yet have mesh size or number of bins. Set them first." << G4endl
129  << "This Method is ignored." << G4endl;
130  return;
131  }
132  if(verboseLevel > 0) G4cout << "G4VScoringMesh::SetPrimitiveScorer() : "
133  << prs->GetName() << " is registered."
134  << " 3D size: ("
135  << fNSegment[0] << ", "
136  << fNSegment[1] << ", "
137  << fNSegment[2] << ")" << G4endl;
138 
139  prs->SetNijk(fNSegment[0], fNSegment[1], fNSegment[2]);
140  fCurrentPS = prs;
141  fMFD->RegisterPrimitive(prs);
143  fMap[prs->GetName()] = map;
144 }
145 
147 
148  if(!fCurrentPS) {
149  G4cerr << "ERROR : G4VScoringMesh::SetSDFilter() : a quantity must be defined first. This method is ignored." << G4endl;
150  return;
151  }
152  if(verboseLevel > 0) G4cout << "G4VScoringMesh::SetFilter() : "
153  << filter->GetName()
154  << " is set to "
155  << fCurrentPS->GetName() << G4endl;
156 
157  G4VSDFilter* oldFilter = fCurrentPS->GetFilter();
158  if(oldFilter)
159  {
160  G4cout << "WARNING : G4VScoringMesh::SetFilter() : " << oldFilter->GetName()
161  << " is overwritten by " << filter->GetName() << G4endl;
162  }
163  fCurrentPS->SetFilter(filter);
164 }
165 
168  if(!fCurrentPS) {
169  G4cerr << "ERROR : G4VScoringMesh::SetCurrentPrimitiveScorer() : The primitive scorer <"
170  << name << "> does not found." << G4endl;
171  }
172 }
173 
175  MeshScoreMap::iterator itr = fMap.find(psname);
176  if(itr == fMap.end()) return false;
177  return true;
178 }
179 
181  MeshScoreMap::iterator itr = fMap.find(psname);
182  if(itr == fMap.end()) {
183  return G4String("");
184  } else {
185  return GetPrimitiveScorer(psname)->GetUnit();
186  }
187 }
188 
190  G4String unit = "";
191  if(!fCurrentPS) {
192  G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
193  msg += " Current primitive scorer is null.";
194  G4cerr << msg << G4endl;
195  }else{
196  unit = fCurrentPS->GetUnit();
197  }
198  return unit;
199 }
200 
202  if(!fCurrentPS) {
203  G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
204  msg += " Current primitive scorer is null.";
205  G4cerr << msg << G4endl;
206  }else{
207  fCurrentPS->SetUnit(unit);
208  }
209 }
210 
212  MeshScoreMap::iterator itr = fMap.find(psname);
213  if(itr == fMap.end()) {
214  return 1.;
215  } else {
216  return GetPrimitiveScorer(psname)->GetUnitValue();
217  }
218 }
219 
220 void G4VScoringMesh::GetDivisionAxisNames(G4String divisionAxisNames[3]) {
221  for(int i = 0; i < 3; i++) divisionAxisNames[i] = fDivisionAxisNames[i];
222 }
223 
225  if(!fMFD) return nullptr;
226 
228  for(G4int i = 0; i < nps; i++) {
230  if(name == prs->GetName()) return prs;
231  }
232 
233  return nullptr;
234 }
235 void G4VScoringMesh::List() const {
236 
237  G4cout << " # of segments: ("
238  << fNSegment[0] << ", "
239  << fNSegment[1] << ", "
240  << fNSegment[2] << ")"
241  << G4endl;
242  G4cout << " displacement: ("
243  << fCenterPosition.x()/cm << ", "
244  << fCenterPosition.y()/cm << ", "
245  << fCenterPosition.z()/cm << ") [cm]"
246  << G4endl;
247  if(fRotationMatrix != 0) {
248  G4cout << " rotation matrix: "
249  << fRotationMatrix->xx() << " "
250  << fRotationMatrix->xy() << " "
251  << fRotationMatrix->xz() << G4endl
252  << " "
253  << fRotationMatrix->yx() << " "
254  << fRotationMatrix->yy() << " "
255  << fRotationMatrix->yz() << G4endl
256  << " "
257  << fRotationMatrix->zx() << " "
258  << fRotationMatrix->zy() << " "
259  << fRotationMatrix->zz() << G4endl;
260  }
261 
262 
263  G4cout << " registered primitve scorers : " << G4endl;
265  G4VPrimitiveScorer * prs;
266  for(int i = 0; i < nps; i++) {
267  prs = fMFD->GetPrimitive(i);
268  G4cout << " " << i << " " << prs->GetName();
269  if(prs->GetFilter() != 0) G4cout << " with " << prs->GetFilter()->GetName();
270  G4cout << G4endl;
271  }
272 }
273 
275  G4cout << "scoring mesh name: " << fWorldName << G4endl;
276  G4cout << "# of G4THitsMap : " << fMap.size() << G4endl;
277  for(auto mp : fMap)
278  {
279  G4cout << "[" << mp.first << "]" << G4endl;
280  mp.second->PrintAllHits();
281  }
282  G4cout << G4endl;
283 
284 }
285 
286 
287 void G4VScoringMesh::DrawMesh(const G4String& psName,G4VScoreColorMap* colorMap,G4int axflg)
288 {
289  fDrawPSName = psName;
290  MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
291  if(fMapItr!=fMap.end()) {
292  fDrawUnit = GetPSUnit(psName);
293  fDrawUnitValue = GetPSUnitValue(psName);
294  Draw(fMapItr->second, colorMap,axflg);
295  } else {
296  G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
297  }
298 }
299 
300 void G4VScoringMesh::DrawMesh(const G4String& psName,G4int idxPlane,G4int iColumn,G4VScoreColorMap* colorMap)
301 {
302  fDrawPSName = psName;
303  MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
304  if(fMapItr!=fMap.end()) {
305  fDrawUnit = GetPSUnit(psName);
306  fDrawUnitValue = GetPSUnitValue(psName);
307  DrawColumn(fMapItr->second,colorMap,idxPlane,iColumn);
308  } else {
309  G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
310  }
311 }
312 
314 {
315  G4String psName = map->GetName();
316  MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
317  *(fMapItr->second) += *map;
318 
319  if(verboseLevel > 9) {
320  G4cout << G4endl;
321  G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
322  G4cout << " PS name : " << psName << G4endl;
323  if(fMapItr == fMap.end()) {
324  G4cout << " " << psName << " was not found." << G4endl;
325  } else {
326  G4cout << " map size : " << map->GetSize() << G4endl;
327  map->PrintAllHits();
328  }
329  G4cout << G4endl;
330  }
331 }
332 
334 {
335  G4String psName = map->GetName();
336  MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
337  *(fMapItr->second) += *map;
338 
339  if(verboseLevel > 9) {
340  G4cout << G4endl;
341  G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
342  G4cout << " PS name : " << psName << G4endl;
343  if(fMapItr == fMap.end()) {
344  G4cout << " " << psName << " was not found." << G4endl;
345  } else {
346  G4cout << " map size : " << map->GetSize() << G4endl;
347  map->PrintAllHits();
348  }
349  G4cout << G4endl;
350  }
351 }
352 
354 {
355  if(fConstructed) {
357  SetupGeometry(fWorldPhys);
359  }
360  if(verboseLevel > 0)
361  G4cout << fWorldPhys->GetName() << " --- All quantities are reset."
362  << G4endl;
363  ResetScore();
364  }
365  else {
366  fConstructed = true;
367  SetupGeometry(fWorldPhys);
368  }
369 }
370 
372 {
373  if(fConstructed) {
377  }
378 
379  if(verboseLevel > 0)
380  G4cout << fWorldPhys->GetName() << " --- All quantities are reset." << G4endl;
381  ResetScore();
382 
383  } else {
384  fConstructed = true;
386  }
387 }
388 
390 {
391  const MeshScoreMap scMap = scMesh->GetScoreMap();
392 
393  MeshScoreMap::const_iterator fMapItr = fMap.begin();
394  MeshScoreMap::const_iterator mapItr = scMap.begin();
395  for(; fMapItr != fMap.end(); fMapItr++) {
396  if(verboseLevel > 9) G4cout << "G4VScoringMesh::Merge()" << fMapItr->first << G4endl;
397  *(fMapItr->second) += *(mapItr->second);
398  mapItr++;
399  }
400 }
401 
G4ThreeVector GetSize() const
const G4String & GetUnit() const
MeshScoreMap GetScoreMap() const
virtual void SetupGeometry(G4VPhysicalVolume *fWorldPhys)=0
const XML_Char * name
Definition: expat.h:151
void Construct(G4VPhysicalVolume *fWorldPhys)
G4String GetCurrentPSUnit()
virtual void Draw(RunScore *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0
MeshScoreMap fMap
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector fCenterPosition
G4String GetName() const
void SetNumberOfSegments(G4int nSegment[3])
G4VPrimitiveScorer * GetPrimitive(G4int id) const
double yx() const
double zz() const
void AddNewDetector(G4VSensitiveDetector *aSD)
Definition: G4SDManager.cc:71
void Merge(const G4VScoringMesh *scMesh)
CLHEP::HepRotation G4RotationMatrix
#define G4endl
Definition: G4ios.hh:61
void Accumulate(G4THitsMap< G4double > *map)
void message(RunManager *runmanager)
Definition: ts_scorers.cc:72
void RotateY(G4double delta)
double yz() const
G4String GetPSUnit(const G4String &psname)
double z() const
void GetDivisionAxisNames(G4String divisionAxisNames[3])
G4bool fGeometryHasBeenDestroyed
void DrawMesh(const G4String &psName, G4VScoreColorMap *colorMap, G4int axflg=111)
double xz() const
void SetFilter(G4VSDFilter *f)
virtual void List() const
void RotateZ(G4double delta)
G4VSDFilter * GetFilter() const
G4String GetName() const
Definition: G4VSDFilter.hh:57
double zx() const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4double fDrawUnitValue
double xx() const
void GetNumberOfSegments(G4int nSegment[3])
virtual void PrintAllHits()
Definition: G4THitsMap.hh:489
void SetPrimitiveScorer(G4VPrimitiveScorer *ps)
G4double GetPSUnitValue(const G4String &psname)
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
virtual ~G4VScoringMesh()
G4GLOB_DLL std::ostream G4cerr
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
void SetFilter(G4VSDFilter *filter)
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
void SetSize(G4double size[3])
double zy() const
G4RotationMatrix * fRotationMatrix
int G4int
Definition: G4Types.hh:78
void SetCurrentPSUnit(const G4String &unit)
virtual size_t GetSize() const
Definition: G4THitsMap.hh:135
void SetUnit(const G4String &unit)
double yy() const
G4bool FindPrimitiveScorer(const G4String &psname)
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:92
void WorkerConstruct(G4VPhysicalVolume *fWorldPhys)
double xy() const
static constexpr double cm
Definition: G4SIunits.hh:119
G4GLOB_DLL std::ostream G4cout
double x() const
virtual void DrawColumn(RunScore *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)=0
G4double fSize[3]
void SetCenterPosition(G4double centerPosition[3])
G4double GetUnitValue() const
void SetNijk(G4int i, G4int j, G4int k)
double y() const
G4VPrimitiveScorer * fCurrentPS
void SetCurrentPrimitiveScorer(const G4String &name)
G4LogicalVolume * fMeshElementLogical
G4String fDivisionAxisNames[3]
G4VPrimitiveScorer * GetPrimitiveScorer(const G4String &name)
void RotateX(G4double delta)
G4MultiFunctionalDetector * fMFD
HepRotation & rotateX(double delta)
Definition: Rotation.cc:66
G4VScoringMesh(const G4String &wName)
G4String fDrawPSName
const G4String & GetName() const
G4bool ReadyForQuantity() const
std::map< G4String, RunScore * > MeshScoreMap
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)