Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
G4VScoringMesh Class Referenceabstract

#include <G4VScoringMesh.hh>

Inheritance diagram for G4VScoringMesh:
G4ScoringBox G4ScoringCylinder

Public Member Functions

 G4VScoringMesh (const G4String &wName)
 
virtual ~G4VScoringMesh ()
 
virtual void Construct (G4VPhysicalVolume *fWorldPhys)=0
 
virtual void WorkerConstruct (G4VPhysicalVolume *fWorldPhys)
 
virtual void List () const
 
const G4StringGetWorldName () const
 
G4bool IsActive () const
 
void Activate (G4bool vl=true)
 
MeshShape GetShape () const
 
void Accumulate (G4THitsMap< G4double > *map)
 
void Merge (const G4VScoringMesh *scMesh)
 
void Dump ()
 
void DrawMesh (const G4String &psName, G4VScoreColorMap *colorMap, G4int axflg=111)
 
void DrawMesh (const G4String &psName, G4int idxPlane, G4int iColumn, G4VScoreColorMap *colorMap)
 
virtual void Draw (std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0
 
virtual void DrawColumn (std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)=0
 
void ResetScore ()
 
void SetSize (G4double size[3])
 
G4ThreeVector GetSize () const
 
void SetCenterPosition (G4double centerPosition[3])
 
G4ThreeVector GetTranslation () const
 
void RotateX (G4double delta)
 
void RotateY (G4double delta)
 
void RotateZ (G4double delta)
 
G4RotationMatrix GetRotationMatrix () const
 
void SetNumberOfSegments (G4int nSegment[3])
 
void GetNumberOfSegments (G4int nSegment[3])
 
void SetPrimitiveScorer (G4VPrimitiveScorer *ps)
 
void SetFilter (G4VSDFilter *filter)
 
void SetCurrentPrimitiveScorer (const G4String &name)
 
G4bool FindPrimitiveScorer (const G4String &psname)
 
G4bool IsCurrentPrimitiveScorerNull ()
 
G4String GetPSUnit (const G4String &psname)
 
G4String GetCurrentPSUnit ()
 
void SetCurrentPSUnit (const G4String &unit)
 
G4double GetPSUnitValue (const G4String &psname)
 
void SetDrawPSName (const G4String &psname)
 
void GetDivisionAxisNames (G4String divisionAxisNames[3])
 
void SetNullToCurrentPrimitiveScorer ()
 
void SetVerboseLevel (G4int vl)
 
MeshScoreMap GetScoreMap () const
 
G4bool ReadyForQuantity () const
 
void SetMeshElementLogical (G4LogicalVolume *val)
 
G4LogicalVolumeGetMeshElementLogical () const
 

Protected Member Functions

G4VPrimitiveScorerGetPrimitiveScorer (const G4String &name)
 

Protected Attributes

G4String fWorldName
 
G4VPrimitiveScorerfCurrentPS
 
G4bool fConstructed
 
G4bool fActive
 
MeshShape fShape
 
G4double fSize [3]
 
G4ThreeVector fCenterPosition
 
G4RotationMatrixfRotationMatrix
 
G4int fNSegment [3]
 
std::map< G4String, G4THitsMap
< G4double > * > 
fMap
 
G4MultiFunctionalDetectorfMFD
 
G4int verboseLevel
 
G4bool sizeIsSet
 
G4bool nMeshIsSet
 
G4String fDrawUnit
 
G4double fDrawUnitValue
 
G4String fDrawPSName
 
G4String fDivisionAxisNames [3]
 
G4LogicalVolumefMeshElementLogical
 

Detailed Description

Definition at line 53 of file G4VScoringMesh.hh.

Constructor & Destructor Documentation

G4VScoringMesh::G4VScoringMesh ( const G4String wName)

Definition at line 46 of file G4VScoringMesh.cc.

References G4SDManager::AddNewDetector(), fDivisionAxisNames, fMFD, fNSegment, fSize, and G4SDManager::GetSDMpointer().

47  : fWorldName(wName),fCurrentPS(0),fConstructed(false),fActive(true),
49  verboseLevel(0),sizeIsSet(false),nMeshIsSet(false),
51 {
53 
54  fSize[0] = fSize[1] = fSize[2] = 0.;
55  fNSegment[0] = fNSegment[1] = fNSegment[2] = 1;
57 }
G4MultiFunctionalDetector * fMFD
G4double fSize[3]
G4double fDrawUnitValue
G4RotationMatrix * fRotationMatrix
G4LogicalVolume * fMeshElementLogical
void AddNewDetector(G4VSensitiveDetector *aSD)
Definition: G4SDManager.cc:67
G4VPrimitiveScorer * fCurrentPS
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
G4String fDivisionAxisNames[3]
G4VScoringMesh::~G4VScoringMesh ( )
virtual

Definition at line 59 of file G4VScoringMesh.cc.

60 {
61  ;
62 }

Member Function Documentation

void G4VScoringMesh::Accumulate ( G4THitsMap< G4double > *  map)

Definition at line 314 of file G4VScoringMesh.cc.

References fMap, G4cout, G4endl, G4VHitsCollection::GetName(), G4THitsMap< T >::GetSize(), G4THitsMap< T >::PrintAllHits(), and verboseLevel.

Referenced by G4ScoringManager::Accumulate().

315 {
316  G4String psName = map->GetName();
317  std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
318  *(fMapItr->second) += *map;
319 
320  if(verboseLevel > 9) {
321  G4cout << G4endl;
322  G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
323  G4cout << " PS name : " << psName << G4endl;
324  if(fMapItr == fMap.end()) {
325  G4cout << " "
326  << psName << " was not found." << G4endl;
327  } else {
328  G4cout << " map size : " << map->GetSize() << G4endl;
329  map->PrintAllHits();
330  }
331  G4cout << G4endl;
332  }
333 }
G4GLOB_DLL std::ostream G4cout
virtual size_t GetSize() const
Definition: G4THitsMap.hh:86
virtual void PrintAllHits()
Definition: G4THitsMap.hh:193
#define G4endl
Definition: G4ios.hh:61
std::map< G4String, G4THitsMap< G4double > * > fMap
void G4VScoringMesh::Activate ( G4bool  vl = true)
inline

Definition at line 76 of file G4VScoringMesh.hh.

References fActive.

77  { fActive = vl; }
virtual void G4VScoringMesh::Construct ( G4VPhysicalVolume fWorldPhys)
pure virtual
virtual void G4VScoringMesh::Draw ( std::map< G4int, G4double * > *  map,
G4VScoreColorMap colorMap,
G4int  axflg = 111 
)
pure virtual

Implemented in G4ScoringBox, and G4ScoringCylinder.

Referenced by DrawMesh().

virtual void G4VScoringMesh::DrawColumn ( std::map< G4int, G4double * > *  map,
G4VScoreColorMap colorMap,
G4int  idxProj,
G4int  idxColumn 
)
pure virtual

Implemented in G4ScoringBox, and G4ScoringCylinder.

Referenced by DrawMesh().

void G4VScoringMesh::DrawMesh ( const G4String psName,
G4VScoreColorMap colorMap,
G4int  axflg = 111 
)

Definition at line 288 of file G4VScoringMesh.cc.

References Draw(), fDrawPSName, fDrawUnit, fDrawUnitValue, fMap, G4cerr, G4endl, GetPSUnit(), and GetPSUnitValue().

Referenced by G4VSceneHandler::AddCompound(), and G4ScoringManager::DrawMesh().

289 {
290  fDrawPSName = psName;
291  std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
292  if(fMapItr!=fMap.end()) {
293  fDrawUnit = GetPSUnit(psName);
294  fDrawUnitValue = GetPSUnitValue(psName);
295  Draw(fMapItr->second->GetMap(), colorMap,axflg);
296  } else {
297  G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
298  }
299 }
G4double GetPSUnitValue(const G4String &psname)
G4String GetPSUnit(const G4String &psname)
G4double fDrawUnitValue
virtual void Draw(std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0
#define G4endl
Definition: G4ios.hh:61
std::map< G4String, G4THitsMap< G4double > * > fMap
G4String fDrawPSName
G4GLOB_DLL std::ostream G4cerr
void G4VScoringMesh::DrawMesh ( const G4String psName,
G4int  idxPlane,
G4int  iColumn,
G4VScoreColorMap colorMap 
)

Definition at line 301 of file G4VScoringMesh.cc.

References DrawColumn(), fDrawPSName, fDrawUnit, fDrawUnitValue, fMap, G4cerr, G4endl, GetPSUnit(), and GetPSUnitValue().

302 {
303  fDrawPSName = psName;
304  std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
305  if(fMapItr!=fMap.end()) {
306  fDrawUnit = GetPSUnit(psName);
307  fDrawUnitValue = GetPSUnitValue(psName);
308  DrawColumn(fMapItr->second->GetMap(),colorMap,idxPlane,iColumn);
309  } else {
310  G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
311  }
312 }
G4double GetPSUnitValue(const G4String &psname)
G4String GetPSUnit(const G4String &psname)
G4double fDrawUnitValue
virtual void DrawColumn(std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)=0
#define G4endl
Definition: G4ios.hh:61
std::map< G4String, G4THitsMap< G4double > * > fMap
G4String fDrawPSName
G4GLOB_DLL std::ostream G4cerr
void G4VScoringMesh::Dump ( )

Definition at line 274 of file G4VScoringMesh.cc.

References fMap, fWorldName, G4cout, and G4endl.

274  {
275  G4cout << "scoring mesh name: " << fWorldName << G4endl;
276 
277  G4cout << "# of G4THitsMap : " << fMap.size() << G4endl;
278  std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.begin();
279  for(; itr != fMap.end(); itr++) {
280  G4cout << "[" << itr->first << "]" << G4endl;
281  itr->second->PrintAllHits();
282  }
283  G4cout << G4endl;
284 
285 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
std::map< G4String, G4THitsMap< G4double > * > fMap
G4bool G4VScoringMesh::FindPrimitiveScorer ( const G4String psname)

Definition at line 172 of file G4VScoringMesh.cc.

References fMap.

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

172  {
173  std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.find(psname);;
174  if(itr == fMap.end()) return false;
175  return true;
176 }
std::map< G4String, G4THitsMap< G4double > * > fMap
G4String G4VScoringMesh::GetCurrentPSUnit ( )

Definition at line 187 of file G4VScoringMesh.cc.

References fCurrentPS, G4cerr, G4endl, and G4VPrimitiveScorer::GetUnit().

Referenced by G4ScoreQuantityMessenger::SetNewValue().

187  {
188  G4String unit = "";
189  if(fCurrentPS == 0) {
190  G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
191  msg += " Current primitive scorer is null.";
192  G4cerr << msg << G4endl;
193  }else{
194  unit = fCurrentPS->GetUnit();
195  }
196  return unit;
197 }
G4VPrimitiveScorer * fCurrentPS
#define G4endl
Definition: G4ios.hh:61
const G4String & GetUnit() const
G4GLOB_DLL std::ostream G4cerr
void G4VScoringMesh::GetDivisionAxisNames ( G4String  divisionAxisNames[3])

Definition at line 218 of file G4VScoringMesh.cc.

References fDivisionAxisNames.

Referenced by G4VScoreWriter::DumpAllQuantitiesToFile(), and G4VScoreWriter::DumpQuantityToFile().

218  {
219  for(int i = 0; i < 3; i++) divisionAxisNames[i] = fDivisionAxisNames[i];
220 }
G4String fDivisionAxisNames[3]
G4LogicalVolume* G4VScoringMesh::GetMeshElementLogical ( ) const
inline

Definition at line 197 of file G4VScoringMesh.hh.

References fMeshElementLogical.

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

198  { return fMeshElementLogical; }
G4LogicalVolume * fMeshElementLogical
void G4VScoringMesh::GetNumberOfSegments ( G4int  nSegment[3])

Definition at line 103 of file G4VScoringMesh.cc.

References fNSegment.

Referenced by G4GMocrenFileSceneHandler::AddSolid(), G4ScoreQuantityMessenger::SetNewValue(), and G4VScoreWriter::SetScoringMesh().

103  {
104  for(int i = 0; i < 3; i++) nSegment[i] = fNSegment[i];
105 }
G4VPrimitiveScorer * G4VScoringMesh::GetPrimitiveScorer ( const G4String name)
protected

Definition at line 222 of file G4VScoringMesh.cc.

References fMFD, G4VPrimitiveScorer::GetName(), G4MultiFunctionalDetector::GetNumberOfPrimitives(), and G4MultiFunctionalDetector::GetPrimitive().

Referenced by GetPSUnit(), GetPSUnitValue(), and SetCurrentPrimitiveScorer().

222  {
223  if(fMFD == 0) return 0;
224 
226  for(G4int i = 0; i < nps; i++) {
228  if(name == ps->GetName()) return ps;
229  }
230 
231  return 0;
232 }
G4String GetName() const
G4MultiFunctionalDetector * fMFD
int G4int
Definition: G4Types.hh:78
G4VPrimitiveScorer * GetPrimitive(G4int id) const
G4String G4VScoringMesh::GetPSUnit ( const G4String psname)

Definition at line 178 of file G4VScoringMesh.cc.

References fMap, GetPrimitiveScorer(), and G4VPrimitiveScorer::GetUnit().

Referenced by DrawMesh(), G4VScoreWriter::DumpAllQuantitiesToFile(), and G4VScoreWriter::DumpQuantityToFile().

178  {
179  std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.find(psname);;
180  if(itr == fMap.end()) {
181  return G4String("");
182  } else {
183  return GetPrimitiveScorer(psname)->GetUnit();
184  }
185 }
G4VPrimitiveScorer * GetPrimitiveScorer(const G4String &name)
const G4String & GetUnit() const
std::map< G4String, G4THitsMap< G4double > * > fMap
G4double G4VScoringMesh::GetPSUnitValue ( const G4String psname)

Definition at line 209 of file G4VScoringMesh.cc.

References fMap, GetPrimitiveScorer(), and G4VPrimitiveScorer::GetUnitValue().

Referenced by DrawMesh(), G4VScoreWriter::DumpAllQuantitiesToFile(), and G4VScoreWriter::DumpQuantityToFile().

209  {
210  std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.find(psname);;
211  if(itr == fMap.end()) {
212  return 1.;
213  } else {
214  return GetPrimitiveScorer(psname)->GetUnitValue();
215  }
216 }
G4VPrimitiveScorer * GetPrimitiveScorer(const G4String &name)
G4double GetUnitValue() const
std::map< G4String, G4THitsMap< G4double > * > fMap
G4RotationMatrix G4VScoringMesh::GetRotationMatrix ( ) const
inline

Definition at line 114 of file G4VScoringMesh.hh.

References fRotationMatrix, and CLHEP::HepRotation::IDENTITY.

Referenced by G4GMocrenFileSceneHandler::AddSolid().

114  {
115  if(fRotationMatrix) return *fRotationMatrix;
116  else return G4RotationMatrix::IDENTITY;
117  }
static DLL_API const HepRotation IDENTITY
Definition: Rotation.h:369
G4RotationMatrix * fRotationMatrix
MeshScoreMap G4VScoringMesh::GetScoreMap ( ) const
inline
MeshShape G4VScoringMesh::GetShape ( ) const
inline

Definition at line 79 of file G4VScoringMesh.hh.

References fShape.

Referenced by G4ScoreQuantityMessenger::SetNewValue(), and G4ScoringMessenger::SetNewValue().

80  { return fShape; }
G4ThreeVector G4VScoringMesh::GetSize ( ) const

Definition at line 83 of file G4VScoringMesh.cc.

References fSize, and sizeIsSet.

Referenced by G4GMocrenFileSceneHandler::AddSolid(), and G4ScoreQuantityMessenger::SetNewValue().

83  {
84  if(sizeIsSet)
85  return G4ThreeVector(fSize[0], fSize[1], fSize[2]);
86  else
87  return G4ThreeVector(0., 0., 0.);
88 }
CLHEP::Hep3Vector G4ThreeVector
G4double fSize[3]
G4ThreeVector G4VScoringMesh::GetTranslation ( ) const
inline

Definition at line 106 of file G4VScoringMesh.hh.

References fCenterPosition.

Referenced by G4GMocrenFileSceneHandler::AddSolid().

106 {return fCenterPosition;}
G4ThreeVector fCenterPosition
const G4String& G4VScoringMesh::GetWorldName ( ) const
inline
G4bool G4VScoringMesh::IsActive ( ) const
inline

Definition at line 73 of file G4VScoringMesh.hh.

References fActive.

Referenced by G4VSceneHandler::AddCompound(), and G4PSHitsModel::DescribeYourselfTo().

74  { return fActive; }
G4bool G4VScoringMesh::IsCurrentPrimitiveScorerNull ( )
inline

Definition at line 132 of file G4VScoringMesh.hh.

References fCurrentPS.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

132  {
133  if(fCurrentPS == NULL) return true;
134  else return false;
135  }
G4VPrimitiveScorer * fCurrentPS
void G4VScoringMesh::List ( ) const
virtual

Reimplemented in G4ScoringBox, and G4ScoringCylinder.

Definition at line 233 of file G4VScoringMesh.cc.

References python.hepunit::cm, fCenterPosition, fMFD, fNSegment, fRotationMatrix, G4cout, G4endl, G4VPrimitiveScorer::GetFilter(), G4VSDFilter::GetName(), G4VPrimitiveScorer::GetName(), G4MultiFunctionalDetector::GetNumberOfPrimitives(), G4MultiFunctionalDetector::GetPrimitive(), CLHEP::Hep3Vector::x(), CLHEP::HepRotation::xx(), CLHEP::HepRotation::xy(), CLHEP::HepRotation::xz(), CLHEP::Hep3Vector::y(), CLHEP::HepRotation::yx(), CLHEP::HepRotation::yy(), CLHEP::HepRotation::yz(), CLHEP::Hep3Vector::z(), CLHEP::HepRotation::zx(), CLHEP::HepRotation::zy(), and CLHEP::HepRotation::zz().

Referenced by G4ScoringCylinder::List(), and G4ScoringBox::List().

233  {
234 
235  G4cout << " # of segments: ("
236  << fNSegment[0] << ", "
237  << fNSegment[1] << ", "
238  << fNSegment[2] << ")"
239  << G4endl;
240  G4cout << " displacement: ("
241  << fCenterPosition.x()/cm << ", "
242  << fCenterPosition.y()/cm << ", "
243  << fCenterPosition.z()/cm << ") [cm]"
244  << G4endl;
245  if(fRotationMatrix != 0) {
246  G4cout << " rotation matrix: "
247  << fRotationMatrix->xx() << " "
248  << fRotationMatrix->xy() << " "
249  << fRotationMatrix->xz() << G4endl
250  << " "
251  << fRotationMatrix->yx() << " "
252  << fRotationMatrix->yy() << " "
253  << fRotationMatrix->yz() << G4endl
254  << " "
255  << fRotationMatrix->zx() << " "
256  << fRotationMatrix->zy() << " "
257  << fRotationMatrix->zz() << G4endl;
258  }
259 
260 
261  G4cout << " registered primitve scorers : " << G4endl;
263  G4VPrimitiveScorer * ps;
264  for(int i = 0; i < nps; i++) {
265  ps = fMFD->GetPrimitive(i);
266  G4cout << " " << i << " " << ps->GetName();
267  if(ps->GetFilter() != 0) G4cout << " with " << ps->GetFilter()->GetName();
268  G4cout << G4endl;
269  }
270 
271 
272 }
double xx() const
G4ThreeVector fCenterPosition
G4String GetName() const
double x() const
G4MultiFunctionalDetector * fMFD
double yy() const
double xz() const
G4String GetName() const
Definition: G4VSDFilter.hh:57
int G4int
Definition: G4Types.hh:78
double zx() const
double z() const
double yz() const
G4GLOB_DLL std::ostream G4cout
G4VPrimitiveScorer * GetPrimitive(G4int id) const
G4RotationMatrix * fRotationMatrix
double zy() const
double y() const
#define G4endl
Definition: G4ios.hh:61
double yx() const
double zz() const
double xy() const
G4VSDFilter * GetFilter() const
void G4VScoringMesh::Merge ( const G4VScoringMesh scMesh)

Definition at line 349 of file G4VScoringMesh.cc.

References fMap, G4cout, G4endl, GetScoreMap(), and verboseLevel.

Referenced by G4ScoringManager::Merge().

350 {
351  const MeshScoreMap scMap = scMesh->GetScoreMap();
352 
353  MeshScoreMap::const_iterator fMapItr = fMap.begin();
354  MeshScoreMap::const_iterator mapItr = scMap.begin();
355  for(; fMapItr != fMap.end(); fMapItr++) {
356  if(verboseLevel > 9) G4cout << "G4VScoringMesh::Merge()" << fMapItr->first << G4endl;
357  *(fMapItr->second) += *(mapItr->second);
358  mapItr++;
359  }
360 }
std::map< G4String, G4THitsMap< G4double > * > MeshScoreMap
G4GLOB_DLL std::ostream G4cout
MeshScoreMap GetScoreMap() const
#define G4endl
Definition: G4ios.hh:61
std::map< G4String, G4THitsMap< G4double > * > fMap
G4bool G4VScoringMesh::ReadyForQuantity ( ) const
inline

Definition at line 159 of file G4VScoringMesh.hh.

References nMeshIsSet, and sizeIsSet.

Referenced by SetPrimitiveScorer().

160  { return (sizeIsSet && nMeshIsSet); }
void G4VScoringMesh::ResetScore ( )

Definition at line 64 of file G4VScoringMesh.cc.

References fMap, G4cout, G4endl, and verboseLevel.

Referenced by G4ScoringCylinder::Construct(), G4ScoringBox::Construct(), and WorkerConstruct().

64  {
65  if(verboseLevel > 9) G4cout << "G4VScoringMesh::ResetScore() is called." << G4endl;
66  std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.begin();
67  for(; itr != fMap.end(); itr++) {
68  if(verboseLevel > 9) G4cout << "G4VScoringMesh::ResetScore()" << itr->first << G4endl;
69  itr->second->clear();
70  }
71 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
std::map< G4String, G4THitsMap< G4double > * > fMap
void G4VScoringMesh::RotateX ( G4double  delta)

Definition at line 106 of file G4VScoringMesh.cc.

References fRotationMatrix, and CLHEP::HepRotation::rotateX().

Referenced by G4ScoringMessenger::SetNewValue().

106  {
108  fRotationMatrix->rotateX(delta);
109 }
HepRotation & rotateX(double delta)
Definition: Rotation.cc:66
CLHEP::HepRotation G4RotationMatrix
G4RotationMatrix * fRotationMatrix
void G4VScoringMesh::RotateY ( G4double  delta)

Definition at line 111 of file G4VScoringMesh.cc.

References fRotationMatrix, and CLHEP::HepRotation::rotateY().

Referenced by G4ScoringMessenger::SetNewValue().

111  {
113  fRotationMatrix->rotateY(delta);
114 }
CLHEP::HepRotation G4RotationMatrix
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
G4RotationMatrix * fRotationMatrix
void G4VScoringMesh::RotateZ ( G4double  delta)

Definition at line 116 of file G4VScoringMesh.cc.

References fRotationMatrix, and CLHEP::HepRotation::rotateZ().

Referenced by G4ScoringMessenger::SetNewValue().

116  {
118  fRotationMatrix->rotateZ(delta);
119 }
CLHEP::HepRotation G4RotationMatrix
G4RotationMatrix * fRotationMatrix
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:92
void G4VScoringMesh::SetCenterPosition ( G4double  centerPosition[3])

Definition at line 89 of file G4VScoringMesh.cc.

References fCenterPosition.

Referenced by G4ScoringMessenger::SetNewValue().

89  {
90  fCenterPosition = G4ThreeVector(centerPosition[0], centerPosition[1], centerPosition[2]);
91 }
G4ThreeVector fCenterPosition
CLHEP::Hep3Vector G4ThreeVector
void G4VScoringMesh::SetCurrentPrimitiveScorer ( const G4String name)

Definition at line 164 of file G4VScoringMesh.cc.

References fCurrentPS, G4cerr, G4endl, and GetPrimitiveScorer().

Referenced by G4ScoreQuantityMessenger::SetNewValue().

164  {
166  if(fCurrentPS == 0) {
167  G4cerr << "ERROR : G4VScoringMesh::SetCurrentPrimitiveScorer() : The primitive scorer <"
168  << name << "> does not found." << G4endl;
169  }
170 }
G4VPrimitiveScorer * GetPrimitiveScorer(const G4String &name)
G4VPrimitiveScorer * fCurrentPS
#define G4endl
Definition: G4ios.hh:61
G4GLOB_DLL std::ostream G4cerr
void G4VScoringMesh::SetCurrentPSUnit ( const G4String unit)

Definition at line 199 of file G4VScoringMesh.cc.

References fCurrentPS, G4cerr, G4endl, and G4VPrimitiveScorer::SetUnit().

Referenced by G4ScoreQuantityMessenger::SetNewValue().

199  {
200  if(fCurrentPS == 0) {
201  G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
202  msg += " Current primitive scorer is null.";
203  G4cerr << msg << G4endl;
204  }else{
205  fCurrentPS->SetUnit(unit);
206  }
207 }
G4VPrimitiveScorer * fCurrentPS
#define G4endl
Definition: G4ios.hh:61
void SetUnit(const G4String &unit)
G4GLOB_DLL std::ostream G4cerr
void G4VScoringMesh::SetDrawPSName ( const G4String psname)
inline

Definition at line 145 of file G4VScoringMesh.hh.

References fDrawPSName.

145 {fDrawPSName = psname;}
G4String fDrawPSName
void G4VScoringMesh::SetFilter ( G4VSDFilter filter)

Definition at line 144 of file G4VScoringMesh.cc.

References fCurrentPS, G4cerr, G4cout, G4endl, G4VPrimitiveScorer::GetFilter(), G4VSDFilter::GetName(), G4VPrimitiveScorer::GetName(), G4VPrimitiveScorer::SetFilter(), and verboseLevel.

Referenced by G4ScoreQuantityMessenger::FParticleCommand(), G4ScoreQuantityMessenger::FParticleWithEnergyCommand(), and G4ScoreQuantityMessenger::SetNewValue().

144  {
145 
146  if(fCurrentPS == 0) {
147  G4cerr << "ERROR : G4VScoringMesh::SetSDFilter() : a quantity must be defined first. This method is ignored." << G4endl;
148  return;
149  }
150  if(verboseLevel > 0) G4cout << "G4VScoringMesh::SetFilter() : "
151  << filter->GetName()
152  << " is set to "
153  << fCurrentPS->GetName() << G4endl;
154 
155  G4VSDFilter* oldFilter = fCurrentPS->GetFilter();
156  if(oldFilter)
157  {
158  G4cout << "WARNING : G4VScoringMesh::SetFilter() : " << oldFilter->GetName()
159  << " is overwritten by " << filter->GetName() << G4endl;
160  }
161  fCurrentPS->SetFilter(filter);
162 }
G4String GetName() const
void SetFilter(G4VSDFilter *f)
G4String GetName() const
Definition: G4VSDFilter.hh:57
G4GLOB_DLL std::ostream G4cout
G4VPrimitiveScorer * fCurrentPS
#define G4endl
Definition: G4ios.hh:61
G4GLOB_DLL std::ostream G4cerr
G4VSDFilter * GetFilter() const
void G4VScoringMesh::SetMeshElementLogical ( G4LogicalVolume val)
inline

Definition at line 195 of file G4VScoringMesh.hh.

References fMeshElementLogical.

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

196  { fMeshElementLogical = val; }
G4LogicalVolume * fMeshElementLogical
void G4VScoringMesh::SetNullToCurrentPrimitiveScorer ( )
inline

Definition at line 151 of file G4VScoringMesh.hh.

References fCurrentPS.

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

151 {fCurrentPS = NULL;}
G4VPrimitiveScorer * fCurrentPS
void G4VScoringMesh::SetNumberOfSegments ( G4int  nSegment[3])

Definition at line 92 of file G4VScoringMesh.cc.

References fNSegment, G4Exception(), JustWarning, and nMeshIsSet.

Referenced by G4ScoringMessenger::MeshBinCommand().

92  {
93  if ( !nMeshIsSet ){
94  for(int i = 0; i < 3; i++) fNSegment[i] = nSegment[i];
95  nMeshIsSet = true;
96  } else {
97  G4String message = " The size of scoring segments can not be changed.";
98  G4Exception("G4VScoringMesh::SetNumberOfSegments()",
99  "DigiHitsUtilsScoreVScoringMesh000", JustWarning,
100  message);
101  }
102 }
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void G4VScoringMesh::SetPrimitiveScorer ( G4VPrimitiveScorer ps)

Definition at line 121 of file G4VScoringMesh.cc.

References fCurrentPS, fMap, fMFD, fNSegment, fWorldName, G4cerr, G4cout, G4endl, G4VPrimitiveScorer::GetName(), ReadyForQuantity(), G4MultiFunctionalDetector::RegisterPrimitive(), G4VPrimitiveScorer::SetNijk(), and verboseLevel.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

121  {
122 
123  if(!ReadyForQuantity())
124  {
125  G4cerr << "ERROR : G4VScoringMesh::SetPrimitiveScorer() : " << ps->GetName()
126  << " does not yet have mesh size or number of bins. Set them first." << G4endl
127  << "This Method is ignored." << G4endl;
128  return;
129  }
130  if(verboseLevel > 0) G4cout << "G4VScoringMesh::SetPrimitiveScorer() : "
131  << ps->GetName() << " is registered."
132  << " 3D size: ("
133  << fNSegment[0] << ", "
134  << fNSegment[1] << ", "
135  << fNSegment[2] << ")" << G4endl;
136 
137  ps->SetNijk(fNSegment[0], fNSegment[1], fNSegment[2]);
138  fCurrentPS = ps;
139  fMFD->RegisterPrimitive(ps);
141  fMap[ps->GetName()] = map;
142 }
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
G4String GetName() const
G4MultiFunctionalDetector * fMFD
void SetNijk(G4int i, G4int j, G4int k)
G4GLOB_DLL std::ostream G4cout
G4bool ReadyForQuantity() const
G4VPrimitiveScorer * fCurrentPS
#define G4endl
Definition: G4ios.hh:61
std::map< G4String, G4THitsMap< G4double > * > fMap
G4GLOB_DLL std::ostream G4cerr
void G4VScoringMesh::SetSize ( G4double  size[3])

Definition at line 72 of file G4VScoringMesh.cc.

References fSize, G4Exception(), JustWarning, and sizeIsSet.

Referenced by G4ScoringMessenger::SetNewValue().

72  {
73  if ( !sizeIsSet ){
74  for(int i = 0; i < 3; i++) fSize[i] = size[i];
75  sizeIsSet = true;
76  }else{
77  G4String message = " The size of scoring mesh can not be changed.";
78  G4Exception("G4VScoringMesh::SetSize()",
79  "DigiHitsUtilsScoreVScoringMesh000", JustWarning,
80  message);
81  }
82 }
G4double fSize[3]
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void G4VScoringMesh::SetVerboseLevel ( G4int  vl)
inline

Definition at line 153 of file G4VScoringMesh.hh.

References verboseLevel.

Referenced by G4ScoringManager::RegisterScoringMesh().

154  { verboseLevel = vl; }
void G4VScoringMesh::WorkerConstruct ( G4VPhysicalVolume fWorldPhys)
virtual

Definition at line 335 of file G4VScoringMesh.cc.

References fConstructed, fMeshElementLogical, fMFD, G4cout, G4endl, G4VPhysicalVolume::GetName(), ResetScore(), G4LogicalVolume::SetSensitiveDetector(), and verboseLevel.

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

336 {
337  if(fConstructed) {
338 
339  if(verboseLevel > 0)
340  G4cout << fWorldPhys->GetName() << " --- All quantities are reset." << G4endl;
341  ResetScore();
342 
343  } else {
344  fConstructed = true;
346  }
347 }
G4MultiFunctionalDetector * fMFD
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
G4LogicalVolume * fMeshElementLogical
#define G4endl
Definition: G4ios.hh:61
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)

Field Documentation

G4bool G4VScoringMesh::fActive
protected

Definition at line 170 of file G4VScoringMesh.hh.

Referenced by Activate(), and IsActive().

G4ThreeVector G4VScoringMesh::fCenterPosition
protected
G4bool G4VScoringMesh::fConstructed
protected
G4VPrimitiveScorer* G4VScoringMesh::fCurrentPS
protected
G4String G4VScoringMesh::fDivisionAxisNames[3]
protected
G4String G4VScoringMesh::fDrawPSName
protected
G4String G4VScoringMesh::fDrawUnit
protected
G4double G4VScoringMesh::fDrawUnitValue
protected
std::map<G4String, G4THitsMap<G4double>* > G4VScoringMesh::fMap
protected
G4LogicalVolume* G4VScoringMesh::fMeshElementLogical
protected
G4MultiFunctionalDetector* G4VScoringMesh::fMFD
protected
G4int G4VScoringMesh::fNSegment[3]
protected
G4RotationMatrix* G4VScoringMesh::fRotationMatrix
protected
MeshShape G4VScoringMesh::fShape
protected
G4double G4VScoringMesh::fSize[3]
protected
G4String G4VScoringMesh::fWorldName
protected
G4bool G4VScoringMesh::nMeshIsSet
protected

Definition at line 184 of file G4VScoringMesh.hh.

Referenced by ReadyForQuantity(), and SetNumberOfSegments().

G4bool G4VScoringMesh::sizeIsSet
protected

Definition at line 183 of file G4VScoringMesh.hh.

Referenced by GetSize(), ReadyForQuantity(), and SetSize().

G4int G4VScoringMesh::verboseLevel
protected

The documentation for this class was generated from the following files: