Geant4 10.7.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 Types

enum class  MeshShape {
  box , cylinder , sphere , realWorldLogVol ,
  probe , undefined = -1
}
 
using EventScore = G4THitsMap< G4double >
 
using RunScore = G4THitsMap< G4StatDouble >
 
using MeshScoreMap = std::map< G4String, RunScore * >
 

Public Member Functions

 G4VScoringMesh (const G4String &wName)
 
virtual ~G4VScoringMesh ()
 
void Construct (G4VPhysicalVolume *fWorldPhys)
 
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 Accumulate (G4THitsMap< G4StatDouble > *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 (RunScore *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0
 
virtual void DrawColumn (RunScore *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
 
G4VPrimitiveScorerGetPrimitiveScorer (const G4String &name)
 
void SetMeshElementLogical (G4LogicalVolume *val)
 
G4LogicalVolumeGetMeshElementLogical () const
 
void SetParallelWorldProcess (G4ParallelWorldProcess *proc)
 
G4ParallelWorldProcessGetParallelWorldProcess () const
 
void GeometryHasBeenDestroyed ()
 
void SetCopyNumberLevel (G4int val)
 
G4int GetCopyNumberLevel () const
 
G4bool LayeredMassFlg ()
 

Protected Member Functions

virtual void SetupGeometry (G4VPhysicalVolume *fWorldPhys)=0
 

Protected Attributes

G4String fWorldName
 
G4VPrimitiveScorerfCurrentPS
 
G4bool fConstructed
 
G4bool fActive
 
MeshShape fShape
 
G4double fSize [3]
 
G4ThreeVector fCenterPosition
 
G4RotationMatrixfRotationMatrix
 
G4int fNSegment [3]
 
MeshScoreMap fMap
 
G4MultiFunctionalDetectorfMFD
 
G4int verboseLevel
 
G4bool sizeIsSet
 
G4bool nMeshIsSet
 
G4String fDrawUnit
 
G4double fDrawUnitValue
 
G4String fDrawPSName
 
G4String fDivisionAxisNames [3]
 
G4LogicalVolumefMeshElementLogical
 
G4ParallelWorldProcessfParallelWorldProcess
 
G4bool fGeometryHasBeenDestroyed
 
G4int copyNumberLevel
 
G4bool layeredMassFlg
 

Detailed Description

Definition at line 53 of file G4VScoringMesh.hh.

Member Typedef Documentation

◆ EventScore

Definition at line 57 of file G4VScoringMesh.hh.

◆ MeshScoreMap

Definition at line 59 of file G4VScoringMesh.hh.

◆ RunScore

Definition at line 58 of file G4VScoringMesh.hh.

Member Enumeration Documentation

◆ MeshShape

enum class G4VScoringMesh::MeshShape
strong
Enumerator
box 
cylinder 
sphere 
realWorldLogVol 
probe 
undefined 

Definition at line 56 of file G4VScoringMesh.hh.

Constructor & Destructor Documentation

◆ G4VScoringMesh()

G4VScoringMesh::G4VScoringMesh ( const G4String wName)

Definition at line 45 of file G4VScoringMesh.cc.

46 : fWorldName(wName),fCurrentPS(nullptr),fConstructed(false),fActive(true),
49 verboseLevel(0),sizeIsSet(false),nMeshIsSet(false),
53{
55
56 fSize[0] = fSize[1] = fSize[2] = 0.;
57 fNSegment[0] = fNSegment[1] = fNSegment[2] = 1;
59}
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:39
void AddNewDetector(G4VSensitiveDetector *aSD)
Definition: G4SDManager.cc:71
G4RotationMatrix * fRotationMatrix
G4double fDrawUnitValue
G4bool fGeometryHasBeenDestroyed
G4MultiFunctionalDetector * fMFD
G4String fDivisionAxisNames[3]
G4LogicalVolume * fMeshElementLogical
G4VPrimitiveScorer * fCurrentPS
G4double fSize[3]
G4ParallelWorldProcess * fParallelWorldProcess

◆ ~G4VScoringMesh()

G4VScoringMesh::~G4VScoringMesh ( )
virtual

Definition at line 61 of file G4VScoringMesh.cc.

61 {
62 ;
63}

Member Function Documentation

◆ Accumulate() [1/2]

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

Definition at line 314 of file G4VScoringMesh.cc.

315{
316 G4String psName = map->GetName();
317 MeshScoreMap::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 << " " << psName << " was not found." << G4endl;
326 } else {
327 G4cout << " map size : " << map->GetSize() << G4endl;
328 map->PrintAllHits();
329 }
330 G4cout << G4endl;
331 }
332}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
MeshScoreMap fMap
virtual size_t GetSize() const
Definition: G4THitsMap.hh:162
virtual void PrintAllHits()
Definition: G4THitsMap.hh:516

◆ Accumulate() [2/2]

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

Definition at line 334 of file G4VScoringMesh.cc.

335{
336 G4String psName = map->GetName();
337 MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
338 *(fMapItr->second) += *map;
339
340 if(verboseLevel > 9) {
341 G4cout << G4endl;
342 G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
343 G4cout << " PS name : " << psName << G4endl;
344 if(fMapItr == fMap.end()) {
345 G4cout << " " << psName << " was not found." << G4endl;
346 } else {
347 G4cout << " map size : " << map->GetSize() << G4endl;
348 map->PrintAllHits();
349 }
350 G4cout << G4endl;
351 }
352}

◆ Activate()

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

Definition at line 88 of file G4VScoringMesh.hh.

89 { fActive = vl; }

◆ Construct()

void G4VScoringMesh::Construct ( G4VPhysicalVolume fWorldPhys)

Definition at line 354 of file G4VScoringMesh.cc.

355{
356 if(fConstructed) {
358 SetupGeometry(fWorldPhys);
360 }
361 if(verboseLevel > 0)
362 G4cout << fWorldName << " --- All quantities are reset."
363 << G4endl;
364 ResetScore();
365 }
366 else {
367 fConstructed = true;
368 SetupGeometry(fWorldPhys);
369 }
370}
virtual void SetupGeometry(G4VPhysicalVolume *fWorldPhys)=0

Referenced by G4RunManager::ConstructScoringWorlds().

◆ Draw()

virtual void G4VScoringMesh::Draw ( RunScore map,
G4VScoreColorMap colorMap,
G4int  axflg = 111 
)
pure virtual

◆ DrawColumn()

virtual void G4VScoringMesh::DrawColumn ( RunScore map,
G4VScoreColorMap colorMap,
G4int  idxProj,
G4int  idxColumn 
)
pure virtual

◆ 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 MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
305 if(fMapItr!=fMap.end()) {
306 fDrawUnit = GetPSUnit(psName);
308 DrawColumn(fMapItr->second,colorMap,idxPlane,iColumn);
309 } else {
310 G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
311 }
312}
G4GLOB_DLL std::ostream G4cerr
G4String fDrawPSName
G4double GetPSUnitValue(const G4String &psname)
G4String GetPSUnit(const G4String &psname)
virtual void DrawColumn(RunScore *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 MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
292 if(fMapItr!=fMap.end()) {
293 fDrawUnit = GetPSUnit(psName);
295 Draw(fMapItr->second, colorMap,axflg);
296 } else {
297 G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
298 }
299}
virtual void Draw(RunScore *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0

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

◆ Dump()

void G4VScoringMesh::Dump ( )

Definition at line 275 of file G4VScoringMesh.cc.

275 {
276 G4cout << "scoring mesh name: " << fWorldName << G4endl;
277 G4cout << "# of G4THitsMap : " << fMap.size() << G4endl;
278 for(auto mp : fMap)
279 {
280 G4cout << "[" << mp.first << "]" << G4endl;
281 mp.second->PrintAllHits();
282 }
283 G4cout << G4endl;
284
285}

◆ FindPrimitiveScorer()

G4bool G4VScoringMesh::FindPrimitiveScorer ( const G4String psname)

Definition at line 175 of file G4VScoringMesh.cc.

175 {
176 MeshScoreMap::iterator itr = fMap.find(psname);
177 if(itr == fMap.end()) return false;
178 return true;
179}

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

◆ GeometryHasBeenDestroyed()

void G4VScoringMesh::GeometryHasBeenDestroyed ( )
inline

◆ GetCopyNumberLevel()

G4int G4VScoringMesh::GetCopyNumberLevel ( ) const
inline

Definition at line 236 of file G4VScoringMesh.hh.

237 { return copyNumberLevel; }

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetCurrentPSUnit()

G4String G4VScoringMesh::GetCurrentPSUnit ( )

Definition at line 190 of file G4VScoringMesh.cc.

190 {
191 G4String unit = "";
192 if(!fCurrentPS) {
193 G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
194 msg += " Current primitive scorer is null.";
195 G4cerr << msg << G4endl;
196 }else{
197 unit = fCurrentPS->GetUnit();
198 }
199 return unit;
200}
const G4String & GetUnit() const

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetDivisionAxisNames()

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

Definition at line 221 of file G4VScoringMesh.cc.

221 {
222 for(int i = 0; i < 3; i++) divisionAxisNames[i] = fDivisionAxisNames[i];
223}

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

◆ GetMeshElementLogical()

G4LogicalVolume * G4VScoringMesh::GetMeshElementLogical ( ) const
inline

Definition at line 212 of file G4VScoringMesh.hh.

213 { return fMeshElementLogical; }

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

◆ GetNumberOfSegments()

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

Definition at line 105 of file G4VScoringMesh.cc.

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

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

◆ GetParallelWorldProcess()

G4ParallelWorldProcess * G4VScoringMesh::GetParallelWorldProcess ( ) const
inline

◆ GetPrimitiveScorer()

G4VPrimitiveScorer * G4VScoringMesh::GetPrimitiveScorer ( const G4String name)

Definition at line 225 of file G4VScoringMesh.cc.

225 {
226 if(!fMFD) return nullptr;
227
229 for(G4int i = 0; i < nps; i++) {
231 if(name == prs->GetName()) return prs;
232 }
233
234 return nullptr;
235}
int G4int
Definition: G4Types.hh:85
G4VPrimitiveScorer * GetPrimitive(G4int id) const
G4String GetName() const

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

◆ GetPSUnit()

G4String G4VScoringMesh::GetPSUnit ( const G4String psname)

Definition at line 181 of file G4VScoringMesh.cc.

181 {
182 MeshScoreMap::iterator itr = fMap.find(psname);
183 if(itr == fMap.end()) {
184 return G4String("");
185 } else {
186 return GetPrimitiveScorer(psname)->GetUnit();
187 }
188}
G4VPrimitiveScorer * GetPrimitiveScorer(const G4String &name)

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

◆ GetPSUnitValue()

G4double G4VScoringMesh::GetPSUnitValue ( const G4String psname)

Definition at line 212 of file G4VScoringMesh.cc.

212 {
213 MeshScoreMap::iterator itr = fMap.find(psname);
214 if(itr == fMap.end()) {
215 return 1.;
216 } else {
217 return GetPrimitiveScorer(psname)->GetUnitValue();
218 }
219}
G4double GetUnitValue() const

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

◆ GetRotationMatrix()

G4RotationMatrix G4VScoringMesh::GetRotationMatrix ( ) const
inline

Definition at line 128 of file G4VScoringMesh.hh.

128 {
130 else return G4RotationMatrix::IDENTITY;
131 }
static DLL_API const HepRotation IDENTITY
Definition: Rotation.h:366

Referenced by G4GMocrenFileSceneHandler::AddSolid().

◆ GetScoreMap()

◆ GetShape()

◆ GetSize()

G4ThreeVector G4VScoringMesh::GetSize ( ) const

Definition at line 85 of file G4VScoringMesh.cc.

85 {
86 if(sizeIsSet)
87 return G4ThreeVector(fSize[0], fSize[1], fSize[2]);
88 else
89 return G4ThreeVector(0., 0., 0.);
90}
CLHEP::Hep3Vector G4ThreeVector

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

◆ GetTranslation()

G4ThreeVector G4VScoringMesh::GetTranslation ( ) const
inline

Definition at line 120 of file G4VScoringMesh.hh.

120{return fCenterPosition;}
G4ThreeVector fCenterPosition

Referenced by G4GMocrenFileSceneHandler::AddSolid().

◆ GetWorldName()

const G4String & G4VScoringMesh::GetWorldName ( ) const
inline

◆ IsActive()

G4bool G4VScoringMesh::IsActive ( ) const
inline

Definition at line 85 of file G4VScoringMesh.hh.

86 { return fActive; }

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

◆ IsCurrentPrimitiveScorerNull()

G4bool G4VScoringMesh::IsCurrentPrimitiveScorerNull ( )
inline

Definition at line 147 of file G4VScoringMesh.hh.

147 {
148 if(fCurrentPS == nullptr) return true;
149 else return false;
150 }

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ LayeredMassFlg()

G4bool G4VScoringMesh::LayeredMassFlg ( )
inline

◆ List()

void G4VScoringMesh::List ( ) const
virtual

Reimplemented in G4ScoringBox, G4ScoringCylinder, G4ScoringProbe, and G4ScoringRealWorld.

Definition at line 236 of file G4VScoringMesh.cc.

236 {
237
238 G4cout << " # of segments: ("
239 << fNSegment[0] << ", "
240 << fNSegment[1] << ", "
241 << fNSegment[2] << ")"
242 << G4endl;
243 G4cout << " displacement: ("
244 << fCenterPosition.x()/cm << ", "
245 << fCenterPosition.y()/cm << ", "
246 << fCenterPosition.z()/cm << ") [cm]"
247 << G4endl;
248 if(fRotationMatrix != 0) {
249 G4cout << " rotation matrix: "
250 << fRotationMatrix->xx() << " "
251 << fRotationMatrix->xy() << " "
252 << fRotationMatrix->xz() << G4endl
253 << " "
254 << fRotationMatrix->yx() << " "
255 << fRotationMatrix->yy() << " "
256 << fRotationMatrix->yz() << G4endl
257 << " "
258 << fRotationMatrix->zx() << " "
259 << fRotationMatrix->zy() << " "
260 << fRotationMatrix->zz() << G4endl;
261 }
262
263
264 G4cout << " registered primitve scorers : " << G4endl;
266 G4VPrimitiveScorer * prs;
267 for(int i = 0; i < nps; i++) {
268 prs = fMFD->GetPrimitive(i);
269 G4cout << " " << i << " " << prs->GetName();
270 if(prs->GetFilter() != 0) G4cout << " with " << prs->GetFilter()->GetName();
271 G4cout << G4endl;
272 }
273}
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:56

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

◆ Merge()

void G4VScoringMesh::Merge ( const G4VScoringMesh scMesh)

Definition at line 390 of file G4VScoringMesh.cc.

391{
392 const MeshScoreMap scMap = scMesh->GetScoreMap();
393
394 MeshScoreMap::const_iterator fMapItr = fMap.begin();
395 MeshScoreMap::const_iterator mapItr = scMap.begin();
396 for(; fMapItr != fMap.end(); fMapItr++) {
397 if(verboseLevel > 9) G4cout << "G4VScoringMesh::Merge()" << fMapItr->first << G4endl;
398 *(fMapItr->second) += *(mapItr->second);
399 mapItr++;
400 }
401}
std::map< G4String, RunScore * > MeshScoreMap
MeshScoreMap GetScoreMap() const

Referenced by G4ScoringManager::Merge().

◆ ReadyForQuantity()

G4bool G4VScoringMesh::ReadyForQuantity ( ) const
inline

Definition at line 174 of file G4VScoringMesh.hh.

175 { return (sizeIsSet && nMeshIsSet); }

Referenced by SetPrimitiveScorer().

◆ ResetScore()

void G4VScoringMesh::ResetScore ( )

Definition at line 65 of file G4VScoringMesh.cc.

65 {
66 if(verboseLevel > 9) G4cout << "G4VScoringMesh::ResetScore() is called." << G4endl;
67 for(auto mp : fMap)
68 {
69 if(verboseLevel > 9) G4cout << "G4VScoringMesh::ResetScore()" << mp.first << G4endl;
70 mp.second->clear();
71 }
72}

Referenced by Construct(), and WorkerConstruct().

◆ RotateX()

void G4VScoringMesh::RotateX ( G4double  delta)

Definition at line 108 of file G4VScoringMesh.cc.

108 {
110 fRotationMatrix->rotateX(delta);
111}
CLHEP::HepRotation G4RotationMatrix
HepRotation & rotateX(double delta)
Definition: Rotation.cc:61

Referenced by G4ScoringMessenger::SetNewValue().

◆ RotateY()

void G4VScoringMesh::RotateY ( G4double  delta)

Definition at line 113 of file G4VScoringMesh.cc.

113 {
115 fRotationMatrix->rotateY(delta);
116}
HepRotation & rotateY(double delta)
Definition: Rotation.cc:74

Referenced by G4ScoringMessenger::SetNewValue().

◆ RotateZ()

void G4VScoringMesh::RotateZ ( G4double  delta)

Definition at line 118 of file G4VScoringMesh.cc.

118 {
120 fRotationMatrix->rotateZ(delta);
121}
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:87

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetCenterPosition()

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

Definition at line 91 of file G4VScoringMesh.cc.

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

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetCopyNumberLevel()

void G4VScoringMesh::SetCopyNumberLevel ( G4int  val)
inline

Definition at line 234 of file G4VScoringMesh.hh.

235 { copyNumberLevel = val; }

◆ SetCurrentPrimitiveScorer()

void G4VScoringMesh::SetCurrentPrimitiveScorer ( const G4String name)

Definition at line 167 of file G4VScoringMesh.cc.

167 {
169 if(!fCurrentPS) {
170 G4cerr << "ERROR : G4VScoringMesh::SetCurrentPrimitiveScorer() : The primitive scorer <"
171 << name << "> does not found." << G4endl;
172 }
173}
const char * name(G4int ptype)

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetCurrentPSUnit()

void G4VScoringMesh::SetCurrentPSUnit ( const G4String unit)

Definition at line 202 of file G4VScoringMesh.cc.

202 {
203 if(!fCurrentPS) {
204 G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
205 msg += " Current primitive scorer is null.";
206 G4cerr << msg << G4endl;
207 }else{
208 fCurrentPS->SetUnit(unit);
209 }
210}
void SetUnit(const G4String &unit)

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetDrawPSName()

void G4VScoringMesh::SetDrawPSName ( const G4String psname)
inline

Definition at line 160 of file G4VScoringMesh.hh.

160{fDrawPSName = psname;}

◆ SetFilter()

void G4VScoringMesh::SetFilter ( G4VSDFilter filter)

Definition at line 147 of file G4VScoringMesh.cc.

147 {
148
149 if(!fCurrentPS) {
150 G4cerr << "ERROR : G4VScoringMesh::SetSDFilter() : a quantity must be defined first. This method is ignored." << G4endl;
151 return;
152 }
153 if(verboseLevel > 0) G4cout << "G4VScoringMesh::SetFilter() : "
154 << filter->GetName()
155 << " is set to "
156 << fCurrentPS->GetName() << G4endl;
157
158 G4VSDFilter* oldFilter = fCurrentPS->GetFilter();
159 if(oldFilter)
160 {
161 G4cout << "WARNING : G4VScoringMesh::SetFilter() : " << oldFilter->GetName()
162 << " is overwritten by " << filter->GetName() << G4endl;
163 }
164 fCurrentPS->SetFilter(filter);
165}
void SetFilter(G4VSDFilter *f)

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

◆ SetMeshElementLogical()

void G4VScoringMesh::SetMeshElementLogical ( G4LogicalVolume val)
inline

Definition at line 210 of file G4VScoringMesh.hh.

211 { fMeshElementLogical = val; }

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

◆ SetNullToCurrentPrimitiveScorer()

void G4VScoringMesh::SetNullToCurrentPrimitiveScorer ( )
inline

Definition at line 166 of file G4VScoringMesh.hh.

166{fCurrentPS = nullptr;}

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

◆ SetNumberOfSegments()

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

Definition at line 94 of file G4VScoringMesh.cc.

94 {
96 for(int i = 0; i < 3; i++) fNSegment[i] = nSegment[i];
97 nMeshIsSet = true;
98 } else {
99 G4String message = " The size of scoring segments can not be changed.";
100 G4Exception("G4VScoringMesh::SetNumberOfSegments()",
101 "DigiHitsUtilsScoreVScoringMesh000", JustWarning,
102 message);
103 }
104}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35

Referenced by G4ScoringProbe::G4ScoringProbe(), G4ScoringRealWorld::G4ScoringRealWorld(), G4ScoringProbe::LocateProbe(), G4ScoringMessenger::MeshBinCommand(), and G4ScoringRealWorld::SetupGeometry().

◆ SetParallelWorldProcess()

void G4VScoringMesh::SetParallelWorldProcess ( G4ParallelWorldProcess proc)
inline

◆ SetPrimitiveScorer()

void G4VScoringMesh::SetPrimitiveScorer ( G4VPrimitiveScorer ps)

Definition at line 123 of file G4VScoringMesh.cc.

123 {
124
125 if(!ReadyForQuantity())
126 {
127 G4cerr << "ERROR : G4VScoringMesh::SetPrimitiveScorer() : "
128 << prs->GetName()
129 << " does not yet have mesh size or number of bins. Set them first." << G4endl
130 << "This Method is ignored." << G4endl;
131 return;
132 }
133 if(verboseLevel > 0) G4cout << "G4VScoringMesh::SetPrimitiveScorer() : "
134 << prs->GetName() << " is registered."
135 << " 3D size: ("
136 << fNSegment[0] << ", "
137 << fNSegment[1] << ", "
138 << fNSegment[2] << ")" << G4endl;
139
140 prs->SetNijk(fNSegment[0], fNSegment[1], fNSegment[2]);
141 fCurrentPS = prs;
144 fMap[prs->GetName()] = map;
145}
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
G4bool ReadyForQuantity() const

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetSize()

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

Definition at line 74 of file G4VScoringMesh.cc.

74 {
75 if ( !sizeIsSet ){
76 for(int i = 0; i < 3; i++) fSize[i] = size[i];
77 sizeIsSet = true;
78 }else{
79 G4String message = " The size of scoring mesh is updated.";
80 G4Exception("G4VScoringMesh::SetSize()",
81 "DigiHitsUtilsScoreVScoringMesh000", JustWarning,
82 message);
83 }
84}

Referenced by G4ScoringProbe::G4ScoringProbe(), G4ScoringRealWorld::G4ScoringRealWorld(), and G4ScoringMessenger::SetNewValue().

◆ SetupGeometry()

virtual void G4VScoringMesh::SetupGeometry ( G4VPhysicalVolume fWorldPhys)
protectedpure virtual

◆ SetVerboseLevel()

void G4VScoringMesh::SetVerboseLevel ( G4int  vl)
inline

Definition at line 168 of file G4VScoringMesh.hh.

169 { verboseLevel = vl; }

Referenced by G4ScoringManager::RegisterScoringMesh().

◆ WorkerConstruct()

void G4VScoringMesh::WorkerConstruct ( G4VPhysicalVolume fWorldPhys)

Definition at line 372 of file G4VScoringMesh.cc.

373{
374 if(fConstructed) {
378 }
379
380 if(verboseLevel > 0)
381 G4cout << fWorldPhys->GetName() << " --- All quantities are reset." << G4endl;
382 ResetScore();
383
384 } else {
385 fConstructed = true;
387 }
388}
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
const G4String & GetName() const

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

Member Data Documentation

◆ copyNumberLevel

G4int G4VScoringMesh::copyNumberLevel
protected

Definition at line 230 of file G4VScoringMesh.hh.

Referenced by GetCopyNumberLevel(), and SetCopyNumberLevel().

◆ fActive

G4bool G4VScoringMesh::fActive
protected

Definition at line 185 of file G4VScoringMesh.hh.

Referenced by Activate(), and IsActive().

◆ fCenterPosition

◆ fConstructed

G4bool G4VScoringMesh::fConstructed
protected

Definition at line 184 of file G4VScoringMesh.hh.

Referenced by Construct(), and WorkerConstruct().

◆ fCurrentPS

◆ fDivisionAxisNames

G4String G4VScoringMesh::fDivisionAxisNames[3]
protected

◆ fDrawPSName

◆ fDrawUnit

G4String G4VScoringMesh::fDrawUnit
protected

◆ fDrawUnitValue

G4double G4VScoringMesh::fDrawUnitValue
protected

◆ fGeometryHasBeenDestroyed

G4bool G4VScoringMesh::fGeometryHasBeenDestroyed
protected

Definition at line 217 of file G4VScoringMesh.hh.

Referenced by Construct(), GeometryHasBeenDestroyed(), and WorkerConstruct().

◆ fMap

◆ fMeshElementLogical

◆ fMFD

◆ fNSegment

◆ fParallelWorldProcess

G4ParallelWorldProcess* G4VScoringMesh::fParallelWorldProcess
protected

Definition at line 216 of file G4VScoringMesh.hh.

Referenced by GetParallelWorldProcess(), and SetParallelWorldProcess().

◆ fRotationMatrix

◆ fShape

◆ fSize

◆ fWorldName

◆ layeredMassFlg

G4bool G4VScoringMesh::layeredMassFlg
protected

Definition at line 243 of file G4VScoringMesh.hh.

Referenced by LayeredMassFlg(), and G4ScoringProbe::SetMaterial().

◆ nMeshIsSet

G4bool G4VScoringMesh::nMeshIsSet
protected

Definition at line 199 of file G4VScoringMesh.hh.

Referenced by ReadyForQuantity(), and SetNumberOfSegments().

◆ sizeIsSet

G4bool G4VScoringMesh::sizeIsSet
protected

Definition at line 198 of file G4VScoringMesh.hh.

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

◆ verboseLevel


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