Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VScoringMesh Class Referenceabstract

#include <G4VScoringMesh.hh>

+ Inheritance diagram for G4VScoringMesh:

Public Member Functions

 G4VScoringMesh (const G4String &wName)
 
virtual ~G4VScoringMesh ()
 
virtual void Construct (G4VPhysicalVolume *fWorldPhys)=0
 
virtual void List () const
 
const G4StringGetWorldName () const
 
G4bool IsActive () const
 
void Activate (G4bool vl=true)
 
MeshShape GetShape () const
 
void Accumulate (G4THitsMap< G4double > *map)
 
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 ()
 
G4bool ReadyForQuantity () 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]
 

Detailed Description

Definition at line 52 of file G4VScoringMesh.hh.

Constructor & Destructor Documentation

◆ G4VScoringMesh()

G4VScoringMesh::G4VScoringMesh ( const G4String wName)

Definition at line 46 of file G4VScoringMesh.cc.

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}
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
void AddNewDetector(G4VSensitiveDetector *aSD)
Definition: G4SDManager.cc:67
G4RotationMatrix * fRotationMatrix
G4double fDrawUnitValue
G4MultiFunctionalDetector * fMFD
G4String fDivisionAxisNames[3]
G4VPrimitiveScorer * fCurrentPS
G4double fSize[3]

◆ ~G4VScoringMesh()

G4VScoringMesh::~G4VScoringMesh ( )
virtual

Definition at line 59 of file G4VScoringMesh.cc.

60{
61 ;
62}

Member Function Documentation

◆ Accumulate()

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

Definition at line 186 of file G4VScoringMesh.hh.

187{
188 G4String psName = map->GetName();
189 std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
190 *(fMapItr->second) += *map;
191
192 if(verboseLevel > 9) {
193 G4cout << G4endl;
194 G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
195 G4cout << " PS name : " << psName << G4endl;
196 if(fMapItr == fMap.end()) {
197 G4cout << " "
198 << psName << " was not found." << G4endl;
199 } else {
200 G4cout << " map size : " << map->GetSize() << G4endl;
201 map->PrintAllHits();
202 }
203 G4cout << G4endl;
204 }
205}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
virtual size_t GetSize() const
Definition: G4THitsMap.hh:86
virtual void PrintAllHits()
Definition: G4THitsMap.hh:193
std::map< G4String, G4THitsMap< G4double > * > fMap

◆ Activate()

void G4VScoringMesh::Activate ( G4bool  vl = true)
inline

Definition at line 72 of file G4VScoringMesh.hh.

73 { fActive = vl; }

◆ Construct()

virtual void G4VScoringMesh::Construct ( G4VPhysicalVolume fWorldPhys)
pure virtual

◆ Draw()

virtual void G4VScoringMesh::Draw ( std::map< G4int, G4double * > *  map,
G4VScoreColorMap colorMap,
G4int  axflg = 111 
)
pure virtual

Implemented in G4ScoringBox, and G4ScoringCylinder.

Referenced by DrawMesh().

◆ DrawColumn()

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().

◆ DrawMesh() [1/2]

void G4VScoringMesh::DrawMesh ( const G4String psName,
G4int  idxPlane,
G4int  iColumn,
G4VScoreColorMap colorMap 
)

Definition at line 301 of file G4VScoringMesh.cc.

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);
308 DrawColumn(fMapItr->second->GetMap(),colorMap,idxPlane,iColumn);
309 } else {
310 G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
311 }
312}
G4DLLIMPORT std::ostream G4cerr
G4String fDrawPSName
G4double GetPSUnitValue(const G4String &psname)
G4String GetPSUnit(const G4String &psname)
virtual void DrawColumn(std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)=0

◆ DrawMesh() [2/2]

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

Definition at line 288 of file G4VScoringMesh.cc.

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);
295 Draw(fMapItr->second->GetMap(), colorMap,axflg);
296 } else {
297 G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
298 }
299}
virtual void Draw(std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0

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

◆ Dump()

void G4VScoringMesh::Dump ( )

Definition at line 274 of file G4VScoringMesh.cc.

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}

◆ FindPrimitiveScorer()

G4bool G4VScoringMesh::FindPrimitiveScorer ( const G4String psname)

Definition at line 172 of file G4VScoringMesh.cc.

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

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

◆ GetCurrentPSUnit()

G4String G4VScoringMesh::GetCurrentPSUnit ( )

Definition at line 187 of file G4VScoringMesh.cc.

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}
const G4String & GetUnit() const

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetDivisionAxisNames()

void G4VScoringMesh::GetDivisionAxisNames ( G4String  divisionAxisNames[3])

Definition at line 218 of file G4VScoringMesh.cc.

218 {
219 for(int i = 0; i < 3; i++) divisionAxisNames[i] = fDivisionAxisNames[i];
220}

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

◆ GetNumberOfSegments()

void G4VScoringMesh::GetNumberOfSegments ( G4int  nSegment[3])

Definition at line 103 of file G4VScoringMesh.cc.

103 {
104 for(int i = 0; i < 3; i++) nSegment[i] = fNSegment[i];
105}

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

◆ GetPrimitiveScorer()

G4VPrimitiveScorer * G4VScoringMesh::GetPrimitiveScorer ( const G4String name)
protected

Definition at line 222 of file G4VScoringMesh.cc.

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}
int G4int
Definition: G4Types.hh:66
G4VPrimitiveScorer * GetPrimitive(G4int id) const
G4String GetName() const

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

◆ GetPSUnit()

G4String G4VScoringMesh::GetPSUnit ( const G4String psname)

Definition at line 178 of file G4VScoringMesh.cc.

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)

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

◆ GetPSUnitValue()

G4double G4VScoringMesh::GetPSUnitValue ( const G4String psname)

Definition at line 209 of file G4VScoringMesh.cc.

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}
G4double GetUnitValue() const

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

◆ GetRotationMatrix()

G4RotationMatrix G4VScoringMesh::GetRotationMatrix ( ) const
inline

Definition at line 108 of file G4VScoringMesh.hh.

108 {
110 else return G4RotationMatrix::IDENTITY;
111 }
static DLL_API const HepRotation IDENTITY
Definition: Rotation.h:369

Referenced by G4GMocrenFileSceneHandler::AddSolid().

◆ GetScoreMap()

◆ GetShape()

MeshShape G4VScoringMesh::GetShape ( ) const
inline

Definition at line 75 of file G4VScoringMesh.hh.

76 { return fShape; }

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

◆ GetSize()

G4ThreeVector G4VScoringMesh::GetSize ( ) const

Definition at line 83 of file G4VScoringMesh.cc.

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

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

◆ GetTranslation()

G4ThreeVector G4VScoringMesh::GetTranslation ( ) const
inline

Definition at line 100 of file G4VScoringMesh.hh.

100{return fCenterPosition;}
G4ThreeVector fCenterPosition

Referenced by G4GMocrenFileSceneHandler::AddSolid().

◆ GetWorldName()

const G4String & G4VScoringMesh::GetWorldName ( ) const
inline

◆ IsActive()

G4bool G4VScoringMesh::IsActive ( ) const
inline

Definition at line 69 of file G4VScoringMesh.hh.

70 { return fActive; }

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

◆ IsCurrentPrimitiveScorerNull()

G4bool G4VScoringMesh::IsCurrentPrimitiveScorerNull ( )
inline

Definition at line 126 of file G4VScoringMesh.hh.

126 {
127 if(fCurrentPS == NULL) return true;
128 else return false;
129 }

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ List()

void G4VScoringMesh::List ( ) const
virtual

Reimplemented in G4ScoringBox, and G4ScoringCylinder.

Definition at line 233 of file G4VScoringMesh.cc.

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;
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 z() const
double x() const
double y() const
double zz() const
double yz() const
double zx() const
double yx() const
double zy() const
double xx() const
double yy() const
double xz() const
double xy() const
G4VSDFilter * GetFilter() const
G4String GetName() const
Definition: G4VSDFilter.hh:57

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

◆ ReadyForQuantity()

G4bool G4VScoringMesh::ReadyForQuantity ( ) const
inline

Definition at line 152 of file G4VScoringMesh.hh.

153 { return (sizeIsSet && nMeshIsSet); }

Referenced by SetPrimitiveScorer().

◆ ResetScore()

void G4VScoringMesh::ResetScore ( )

Definition at line 64 of file G4VScoringMesh.cc.

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}

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

◆ RotateX()

void G4VScoringMesh::RotateX ( G4double  delta)

Definition at line 106 of file G4VScoringMesh.cc.

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

Referenced by G4ScoringMessenger::SetNewValue().

◆ RotateY()

void G4VScoringMesh::RotateY ( G4double  delta)

Definition at line 111 of file G4VScoringMesh.cc.

111 {
113 fRotationMatrix->rotateY(delta);
114}
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79

Referenced by G4ScoringMessenger::SetNewValue().

◆ RotateZ()

void G4VScoringMesh::RotateZ ( G4double  delta)

Definition at line 116 of file G4VScoringMesh.cc.

116 {
118 fRotationMatrix->rotateZ(delta);
119}
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:92

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetCenterPosition()

void G4VScoringMesh::SetCenterPosition ( G4double  centerPosition[3])

Definition at line 89 of file G4VScoringMesh.cc.

89 {
90 fCenterPosition = G4ThreeVector(centerPosition[0], centerPosition[1], centerPosition[2]);
91}

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetCurrentPrimitiveScorer()

void G4VScoringMesh::SetCurrentPrimitiveScorer ( const G4String name)

Definition at line 164 of file G4VScoringMesh.cc.

164 {
166 if(fCurrentPS == 0) {
167 G4cerr << "ERROR : G4VScoringMesh::SetCurrentPrimitiveScorer() : The primitive scorer <"
168 << name << "> does not found." << G4endl;
169 }
170}

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetCurrentPSUnit()

void G4VScoringMesh::SetCurrentPSUnit ( const G4String unit)

Definition at line 199 of file G4VScoringMesh.cc.

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}
void SetUnit(const G4String &unit)

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetDrawPSName()

void G4VScoringMesh::SetDrawPSName ( const G4String psname)
inline

Definition at line 139 of file G4VScoringMesh.hh.

139{fDrawPSName = psname;}

◆ SetFilter()

void G4VScoringMesh::SetFilter ( G4VSDFilter filter)

Definition at line 144 of file G4VScoringMesh.cc.

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}
void SetFilter(G4VSDFilter *f)

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

◆ SetNullToCurrentPrimitiveScorer()

void G4VScoringMesh::SetNullToCurrentPrimitiveScorer ( )
inline

Definition at line 145 of file G4VScoringMesh.hh.

145{fCurrentPS = NULL;}

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

◆ SetNumberOfSegments()

void G4VScoringMesh::SetNumberOfSegments ( G4int  nSegment[3])

Definition at line 92 of file G4VScoringMesh.cc.

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}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by G4ScoringMessenger::MeshBinCommand().

◆ SetPrimitiveScorer()

void G4VScoringMesh::SetPrimitiveScorer ( G4VPrimitiveScorer ps)

Definition at line 121 of file G4VScoringMesh.cc.

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;
141 fMap[ps->GetName()] = map;
142}
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
void SetNijk(G4int i, G4int j, G4int k)
G4bool ReadyForQuantity() const

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetSize()

void G4VScoringMesh::SetSize ( G4double  size[3])

Definition at line 72 of file G4VScoringMesh.cc.

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}

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetVerboseLevel()

void G4VScoringMesh::SetVerboseLevel ( G4int  vl)
inline

Definition at line 147 of file G4VScoringMesh.hh.

148 { verboseLevel = vl; }

Referenced by G4ScoringManager::RegisterScoringMesh().

Member Data Documentation

◆ fActive

G4bool G4VScoringMesh::fActive
protected

Definition at line 163 of file G4VScoringMesh.hh.

Referenced by Activate(), and IsActive().

◆ fCenterPosition

◆ fConstructed

G4bool G4VScoringMesh::fConstructed
protected

Definition at line 162 of file G4VScoringMesh.hh.

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

◆ fCurrentPS

◆ fDivisionAxisNames

G4String G4VScoringMesh::fDivisionAxisNames[3]
protected

◆ fDrawPSName

◆ fDrawUnit

G4String G4VScoringMesh::fDrawUnit
protected

◆ fDrawUnitValue

G4double G4VScoringMesh::fDrawUnitValue
protected

◆ fMap

std::map<G4String, G4THitsMap<G4double>* > G4VScoringMesh::fMap
protected

◆ fMFD

G4MultiFunctionalDetector* G4VScoringMesh::fMFD
protected

Definition at line 172 of file G4VScoringMesh.hh.

Referenced by G4VScoringMesh(), GetPrimitiveScorer(), List(), and SetPrimitiveScorer().

◆ fNSegment

◆ fRotationMatrix

◆ fShape

MeshShape G4VScoringMesh::fShape
protected

◆ fSize

◆ fWorldName

G4String G4VScoringMesh::fWorldName
protected

◆ nMeshIsSet

G4bool G4VScoringMesh::nMeshIsSet
protected

Definition at line 177 of file G4VScoringMesh.hh.

Referenced by ReadyForQuantity(), and SetNumberOfSegments().

◆ sizeIsSet

G4bool G4VScoringMesh::sizeIsSet
protected

Definition at line 176 of file G4VScoringMesh.hh.

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

◆ verboseLevel

G4int G4VScoringMesh::verboseLevel
protected

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