Geant4 11.2.2
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 ()=default
 
virtual void Construct (G4VPhysicalVolume *fWorldPhys)
 
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 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 SetAngles (G4double, G4double)
 
G4double GetStartAngle () const
 
G4double GetAngleSpan () 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]
 
G4double fAngle [2]
 
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 54 of file G4VScoringMesh.hh.

Member Typedef Documentation

◆ EventScore

Definition at line 66 of file G4VScoringMesh.hh.

◆ MeshScoreMap

Definition at line 68 of file G4VScoringMesh.hh.

◆ RunScore

Member Enumeration Documentation

◆ MeshShape

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

Definition at line 57 of file G4VScoringMesh.hh.

Constructor & Destructor Documentation

◆ G4VScoringMesh()

G4VScoringMesh::G4VScoringMesh ( const G4String & wName)

Definition at line 45 of file G4VScoringMesh.cc.

46 : fWorldName(wName)
47 , fCurrentPS(nullptr)
48 , fConstructed(false)
49 , fActive(true)
51 , fRotationMatrix(nullptr)
53 , verboseLevel(0)
54 , sizeIsSet(false)
55 , nMeshIsSet(false)
56 , fDrawUnit("")
57 , fDrawUnitValue(1.)
58 , fMeshElementLogical(nullptr)
59 , fParallelWorldProcess(nullptr)
62 , layeredMassFlg(false)
63{
65
66 fSize[0] = fSize[1] = fSize[2] = 0.;
67 fAngle[0] = 0.0;
68 fAngle[1] = CLHEP::twopi * rad;
69 fNSegment[0] = fNSegment[1] = fNSegment[2] = 1;
71}
static G4SDManager * GetSDMpointer()
void AddNewDetector(G4VSensitiveDetector *aSD)
G4RotationMatrix * fRotationMatrix
G4double fDrawUnitValue
G4bool fGeometryHasBeenDestroyed
G4MultiFunctionalDetector * fMFD
G4String fDivisionAxisNames[3]
G4LogicalVolume * fMeshElementLogical
G4VPrimitiveScorer * fCurrentPS
G4double fAngle[2]
G4double fSize[3]
G4ParallelWorldProcess * fParallelWorldProcess

◆ ~G4VScoringMesh()

virtual G4VScoringMesh::~G4VScoringMesh ( )
virtualdefault

Member Function Documentation

◆ Accumulate() [1/2]

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

Definition at line 383 of file G4VScoringMesh.cc.

384{
385 G4String psName = map->GetName();
386 const auto fMapItr = fMap.find(psName);
387 *(fMapItr->second) += *map;
388
389 if(verboseLevel > 9)
390 {
391 G4cout << G4endl;
392 G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
393 G4cout << " PS name : " << psName << G4endl;
394 if(fMapItr == fMap.cend())
395 {
396 G4cout << " " << psName << " was not found." << G4endl;
397 }
398 else
399 {
400 G4cout << " map size : " << map->GetSize() << G4endl;
401 map->PrintAllHits();
402 }
403 G4cout << G4endl;
404 }
405}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
MeshScoreMap fMap

◆ Accumulate() [2/2]

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

Definition at line 407 of file G4VScoringMesh.cc.

408{
409 G4String psName = map->GetName();
410 const auto fMapItr = fMap.find(psName);
411 *(fMapItr->second) += *map;
412
413 if(verboseLevel > 9)
414 {
415 G4cout << G4endl;
416 G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
417 G4cout << " PS name : " << psName << G4endl;
418 if(fMapItr == fMap.cend())
419 {
420 G4cout << " " << psName << " was not found." << G4endl;
421 }
422 else
423 {
424 G4cout << " map size : " << map->GetSize() << G4endl;
425 map->PrintAllHits();
426 }
427 G4cout << G4endl;
428 }
429}

◆ Activate()

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

Definition at line 89 of file G4VScoringMesh.hh.

89{ fActive = vl; }

◆ Construct()

void G4VScoringMesh::Construct ( G4VPhysicalVolume * fWorldPhys)
virtual

Definition at line 431 of file G4VScoringMesh.cc.

432{
433 if(fConstructed)
434 {
436 {
437 SetupGeometry(fWorldPhys);
439 }
440 if(verboseLevel > 0)
441 G4cout << fWorldName << " --- All quantities are reset." << G4endl;
442 ResetScore();
443 }
444 else
445 {
446 fConstructed = true;
447 SetupGeometry(fWorldPhys);
448 }
449}
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 365 of file G4VScoringMesh.cc.

367{
368 fDrawPSName = psName;
369 const auto fMapItr = fMap.find(psName);
370 if(fMapItr != fMap.cend())
371 {
372 fDrawUnit = GetPSUnit(psName);
374 DrawColumn(fMapItr->second, colorMap, idxPlane, iColumn);
375 }
376 else
377 {
378 G4cerr << "Scorer <" << psName << "> is not defined. Method ignored."
379 << G4endl;
380 }
381}
G4GLOB_DLL std::ostream G4cerr
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 347 of file G4VScoringMesh.cc.

349{
350 fDrawPSName = psName;
351 const auto fMapItr = fMap.find(psName);
352 if(fMapItr != fMap.cend())
353 {
354 fDrawUnit = GetPSUnit(psName);
356 Draw(fMapItr->second, colorMap, axflg);
357 }
358 else
359 {
360 G4cerr << "Scorer <" << psName << "> is not defined. Method ignored."
361 << G4endl;
362 }
363}
virtual void Draw(RunScore *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0

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

◆ Dump()

void G4VScoringMesh::Dump ( )

Definition at line 335 of file G4VScoringMesh.cc.

336{
337 G4cout << "scoring mesh name: " << fWorldName << G4endl;
338 G4cout << "# of G4THitsMap : " << fMap.size() << G4endl;
339 for(const auto& mp : fMap)
340 {
341 G4cout << "[" << mp.first << "]" << G4endl;
342 mp.second->PrintAllHits();
343 }
344 G4cout << G4endl;
345}

◆ FindPrimitiveScorer()

G4bool G4VScoringMesh::FindPrimitiveScorer ( const G4String & psname)

Definition at line 225 of file G4VScoringMesh.cc.

226{
227 const auto itr = fMap.find(psname);
228 return itr != fMap.cend();
229}

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

◆ GeometryHasBeenDestroyed()

void G4VScoringMesh::GeometryHasBeenDestroyed ( )
inline

◆ GetAngleSpan()

G4double G4VScoringMesh::GetAngleSpan ( ) const
inline

Definition at line 123 of file G4VScoringMesh.hh.

123{ return fAngle[1]; }

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetCopyNumberLevel()

G4int G4VScoringMesh::GetCopyNumberLevel ( ) const
inline

Definition at line 214 of file G4VScoringMesh.hh.

214{ return copyNumberLevel; }

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetCurrentPSUnit()

G4String G4VScoringMesh::GetCurrentPSUnit ( )

Definition at line 242 of file G4VScoringMesh.cc.

243{
244 G4String unit = "";
245 if(fCurrentPS == nullptr)
246 {
247 G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
248 msg += " Current primitive scorer is null.";
249 G4cerr << msg << G4endl;
250 }
251 else
252 {
253 unit = fCurrentPS->GetUnit();
254 }
255 return unit;
256}
const G4String & GetUnit() const

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetDivisionAxisNames()

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

Definition at line 283 of file G4VScoringMesh.cc.

284{
285 for(G4int i = 0; i < 3; ++i)
286 divisionAxisNames[i] = fDivisionAxisNames[i];
287}
int G4int
Definition G4Types.hh:85

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

◆ GetMeshElementLogical()

G4LogicalVolume * G4VScoringMesh::GetMeshElementLogical ( ) const
inline

Definition at line 192 of file G4VScoringMesh.hh.

193 {
194 return fMeshElementLogical;
195 }

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

◆ GetNumberOfSegments()

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

Definition at line 141 of file G4VScoringMesh.cc.

142{
143 for(G4int i = 0; i < 3; ++i)
144 nSegment[i] = fNSegment[i];
145}

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 289 of file G4VScoringMesh.cc.

290{
291 if(fMFD == nullptr)
292 return nullptr;
293
295 for(G4int i = 0; i < nps; ++i)
296 {
298 if(name == prs->GetName())
299 return prs;
300 }
301
302 return nullptr;
303}
G4VPrimitiveScorer * GetPrimitive(G4int id) const
G4String GetName() const

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

◆ GetPSUnit()

G4String G4VScoringMesh::GetPSUnit ( const G4String & psname)

Definition at line 231 of file G4VScoringMesh.cc.

232{
233 const auto itr = fMap.find(psname);
234 if(itr == fMap.cend())
235 {
236 return G4String("");
237 }
238
239 return GetPrimitiveScorer(psname)->GetUnit();
240}
G4VPrimitiveScorer * GetPrimitiveScorer(const G4String &name)

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

◆ GetPSUnitValue()

G4double G4VScoringMesh::GetPSUnitValue ( const G4String & psname)

Definition at line 272 of file G4VScoringMesh.cc.

273{
274 const auto itr = fMap.find(psname);
275 if(itr == fMap.cend())
276 {
277 return 1.;
278 }
279
280 return GetPrimitiveScorer(psname)->GetUnitValue();
281}
G4double GetUnitValue() const

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

◆ GetRotationMatrix()

G4RotationMatrix G4VScoringMesh::GetRotationMatrix ( ) const
inline

Definition at line 135 of file G4VScoringMesh.hh.

136 {
137 if(fRotationMatrix != nullptr)
138 return *fRotationMatrix;
140 }
static DLL_API const HepRotation IDENTITY
Definition Rotation.h:366

Referenced by G4GMocrenFileSceneHandler::AddSolid().

◆ GetScoreMap()

◆ GetShape()

◆ GetSize()

G4ThreeVector G4VScoringMesh::GetSize ( ) const

Definition at line 104 of file G4VScoringMesh.cc.

105{
106 if(sizeIsSet)
107 return G4ThreeVector(fSize[0], fSize[1], fSize[2]);
108 return G4ThreeVector(0., 0., 0.);
109}
CLHEP::Hep3Vector G4ThreeVector

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

◆ GetStartAngle()

G4double G4VScoringMesh::GetStartAngle ( ) const
inline

Definition at line 122 of file G4VScoringMesh.hh.

122{ return fAngle[0]; }

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetTranslation()

G4ThreeVector G4VScoringMesh::GetTranslation ( ) const
inline

Definition at line 127 of file G4VScoringMesh.hh.

127{ return fCenterPosition; }
G4ThreeVector fCenterPosition

Referenced by G4GMocrenFileSceneHandler::AddSolid().

◆ GetWorldName()

const G4String & G4VScoringMesh::GetWorldName ( ) const
inline

◆ IsActive()

G4bool G4VScoringMesh::IsActive ( ) const
inline

◆ IsCurrentPrimitiveScorerNull()

G4bool G4VScoringMesh::IsCurrentPrimitiveScorerNull ( )
inline

Definition at line 157 of file G4VScoringMesh.hh.

158 {
159 return fCurrentPS == nullptr;
160 }

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 305 of file G4VScoringMesh.cc.

306{
307 G4cout << " # of segments: (" << fNSegment[0] << ", " << fNSegment[1] << ", "
308 << fNSegment[2] << ")" << G4endl;
309 G4cout << " displacement: (" << fCenterPosition.x() / cm << ", "
310 << fCenterPosition.y() / cm << ", " << fCenterPosition.z() / cm
311 << ") [cm]" << G4endl;
312 if(fRotationMatrix != nullptr)
313 {
314 G4cout << " rotation matrix: " << fRotationMatrix->xx() << " "
315 << fRotationMatrix->xy() << " " << fRotationMatrix->xz() << G4endl
316 << " " << fRotationMatrix->yx() << " "
317 << fRotationMatrix->yy() << " " << fRotationMatrix->yz() << G4endl
318 << " " << fRotationMatrix->zx() << " "
319 << fRotationMatrix->zy() << " " << fRotationMatrix->zz() << G4endl;
320 }
321
322 G4cout << " registered primitve scorers : " << G4endl;
325 for(G4int i = 0; i < nps; ++i)
326 {
327 prs = fMFD->GetPrimitive(i);
328 G4cout << " " << i << " " << prs->GetName();
329 if(prs->GetFilter() != nullptr)
330 G4cout << " with " << prs->GetFilter()->GetName();
331 G4cout << G4endl;
332 }
333}
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

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

◆ Merge()

void G4VScoringMesh::Merge ( const G4VScoringMesh * scMesh)

Definition at line 473 of file G4VScoringMesh.cc.

474{
475 const MeshScoreMap scMap = scMesh->GetScoreMap();
476
477 auto fMapItr = fMap.cbegin();
478 auto mapItr = scMap.cbegin();
479 for(; fMapItr != fMap.cend(); ++fMapItr)
480 {
481 if(verboseLevel > 9)
482 G4cout << "G4VScoringMesh::Merge()" << fMapItr->first << G4endl;
483 *(fMapItr->second) += *(mapItr->second);
484 ++mapItr;
485 }
486}
std::map< G4String, RunScore * > MeshScoreMap
MeshScoreMap GetScoreMap() const

Referenced by G4ScoringManager::Merge().

◆ ReadyForQuantity()

G4bool G4VScoringMesh::ReadyForQuantity ( ) const
inline

Definition at line 182 of file G4VScoringMesh.hh.

182{ return (sizeIsSet && nMeshIsSet); }

Referenced by SetPrimitiveScorer().

◆ ResetScore()

void G4VScoringMesh::ResetScore ( )

Definition at line 73 of file G4VScoringMesh.cc.

74{
75 if(verboseLevel > 9)
76 G4cout << "G4VScoringMesh::ResetScore() is called." << G4endl;
77 for(const auto& mp : fMap)
78 {
79 if(verboseLevel > 9)
80 G4cout << "G4VScoringMesh::ResetScore()" << mp.first << G4endl;
81 mp.second->clear();
82 }
83}

Referenced by Construct(), and WorkerConstruct().

◆ RotateX()

void G4VScoringMesh::RotateX ( G4double delta)

Definition at line 147 of file G4VScoringMesh.cc.

148{
149 if(fRotationMatrix == nullptr)
151 fRotationMatrix->rotateX(delta);
152}
CLHEP::HepRotation G4RotationMatrix
HepRotation & rotateX(double delta)
Definition Rotation.cc:61

Referenced by G4ScoringMessenger::SetNewValue().

◆ RotateY()

void G4VScoringMesh::RotateY ( G4double delta)

Definition at line 154 of file G4VScoringMesh.cc.

155{
156 if(fRotationMatrix == nullptr)
158 fRotationMatrix->rotateY(delta);
159}
HepRotation & rotateY(double delta)
Definition Rotation.cc:74

Referenced by G4ScoringMessenger::SetNewValue().

◆ RotateZ()

void G4VScoringMesh::RotateZ ( G4double delta)

Definition at line 161 of file G4VScoringMesh.cc.

162{
163 if(fRotationMatrix == nullptr)
165 fRotationMatrix->rotateZ(delta);
166}
HepRotation & rotateZ(double delta)
Definition Rotation.cc:87

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetAngles()

void G4VScoringMesh::SetAngles ( G4double startAngle,
G4double spanAngle )

Definition at line 111 of file G4VScoringMesh.cc.

112{
113 fAngle[0] = startAngle;
114 fAngle[1] = spanAngle;
115}

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetCenterPosition()

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

Definition at line 117 of file G4VScoringMesh.cc.

118{
120 G4ThreeVector(centerPosition[0], centerPosition[1], centerPosition[2]);
121}

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetCopyNumberLevel()

void G4VScoringMesh::SetCopyNumberLevel ( G4int val)
inline

Definition at line 213 of file G4VScoringMesh.hh.

213{ copyNumberLevel = val; }

◆ SetCurrentPrimitiveScorer()

void G4VScoringMesh::SetCurrentPrimitiveScorer ( const G4String & name)

Definition at line 214 of file G4VScoringMesh.cc.

215{
217 if(fCurrentPS == nullptr)
218 {
219 G4cerr << "ERROR : G4VScoringMesh::SetCurrentPrimitiveScorer() : The "
220 "primitive scorer <"
221 << name << "> does not found." << G4endl;
222 }
223}
const char * name(G4int ptype)

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetCurrentPSUnit()

void G4VScoringMesh::SetCurrentPSUnit ( const G4String & unit)

Definition at line 258 of file G4VScoringMesh.cc.

259{
260 if(fCurrentPS == nullptr)
261 {
262 G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
263 msg += " Current primitive scorer is null.";
264 G4cerr << msg << G4endl;
265 }
266 else
267 {
268 fCurrentPS->SetUnit(unit);
269 }
270}
void SetUnit(const G4String &unit)

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetDrawPSName()

void G4VScoringMesh::SetDrawPSName ( const G4String & psname)
inline

Definition at line 170 of file G4VScoringMesh.hh.

170{ fDrawPSName = psname; }

◆ SetFilter()

void G4VScoringMesh::SetFilter ( G4VSDFilter * filter)

Definition at line 192 of file G4VScoringMesh.cc.

193{
194 if(fCurrentPS == nullptr)
195 {
196 G4cerr << "ERROR : G4VScoringMesh::SetSDFilter() : a quantity must be "
197 "defined first. This method is ignored."
198 << G4endl;
199 return;
200 }
201 if(verboseLevel > 0)
202 G4cout << "G4VScoringMesh::SetFilter() : " << filter->GetName()
203 << " is set to " << fCurrentPS->GetName() << G4endl;
204
205 G4VSDFilter* oldFilter = fCurrentPS->GetFilter();
206 if(oldFilter != nullptr)
207 {
208 G4cout << "WARNING : G4VScoringMesh::SetFilter() : " << oldFilter->GetName()
209 << " is overwritten by " << filter->GetName() << G4endl;
210 }
211 fCurrentPS->SetFilter(filter);
212}
void SetFilter(G4VSDFilter *f)

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

◆ SetMeshElementLogical()

void G4VScoringMesh::SetMeshElementLogical ( G4LogicalVolume * val)
inline

Definition at line 188 of file G4VScoringMesh.hh.

189 {
191 }

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

◆ SetNullToCurrentPrimitiveScorer()

void G4VScoringMesh::SetNullToCurrentPrimitiveScorer ( )
inline

Definition at line 176 of file G4VScoringMesh.hh.

176{ fCurrentPS = nullptr; }

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

◆ SetNumberOfSegments()

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

Definition at line 123 of file G4VScoringMesh.cc.

124{
127 {
128 for(G4int i = 0; i < 3; ++i)
129 fNSegment[i] = nSegment[i];
130 nMeshIsSet = true;
131 }
132 else
133 {
134 G4String message = " Number of bins has already been set and it cannot be changed.\n";
135 message += " This method is ignored.";
136 G4Exception("G4VScoringMesh::SetNumberOfSegments()",
137 "DigiHitsUtilsScoreVScoringMesh000", JustWarning, message);
138 }
139}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)

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 168 of file G4VScoringMesh.cc.

169{
170 if(!ReadyForQuantity())
171 {
172 G4cerr << "ERROR : G4VScoringMesh::SetPrimitiveScorer() : "
173 << prs->GetName()
174 << " does not yet have mesh size or number of bins. Set them first."
175 << G4endl << "This Method is ignored." << G4endl;
176 return;
177 }
178 if(verboseLevel > 0)
179 G4cout << "G4VScoringMesh::SetPrimitiveScorer() : " << prs->GetName()
180 << " is registered."
181 << " 3D size: (" << fNSegment[0] << ", " << fNSegment[1] << ", "
182 << fNSegment[2] << ")" << G4endl;
183
184 prs->SetNijk(fNSegment[0], fNSegment[1], fNSegment[2]);
185 fCurrentPS = prs;
187 auto map =
188 new G4THitsMap<G4StatDouble>(fWorldName, prs->GetName());
189 fMap[prs->GetName()] = map;
190}
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
G4bool ReadyForQuantity() const

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetSize()

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

Definition at line 85 of file G4VScoringMesh.cc.

86{
87 if(!sizeIsSet)
88 {
89 sizeIsSet = true;
90 for(G4int i = 0; i < 3; ++i)
91 {
92 fSize[i] = size[i];
93 }
94 }
95 else
96 {
97 G4String message = " Mesh size has already been set and it cannot be changed.\n";
98 message += " This method is ignored.";
99 G4Exception("G4VScoringMesh::SetSize()",
100 "DigiHitsUtilsScoreVScoringMesh000", JustWarning, message);
101 }
102}

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 178 of file G4VScoringMesh.hh.

178{ verboseLevel = vl; }

Referenced by G4ScoringManager::RegisterScoringMesh().

◆ WorkerConstruct()

void G4VScoringMesh::WorkerConstruct ( G4VPhysicalVolume * fWorldPhys)
virtual

Definition at line 451 of file G4VScoringMesh.cc.

452{
453 if(fConstructed)
454 {
456 {
459 }
460
461 if(verboseLevel > 0)
462 G4cout << fWorldPhys->GetName() << " --- All quantities are reset."
463 << G4endl;
464 ResetScore();
465 }
466 else
467 {
468 fConstructed = true;
470 }
471}
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
const G4String & GetName() const

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

Member Data Documentation

◆ copyNumberLevel

G4int G4VScoringMesh::copyNumberLevel
protected

Definition at line 254 of file G4VScoringMesh.hh.

Referenced by GetCopyNumberLevel(), and SetCopyNumberLevel().

◆ fActive

G4bool G4VScoringMesh::fActive
protected

Definition at line 226 of file G4VScoringMesh.hh.

Referenced by Activate(), and IsActive().

◆ fAngle

◆ fCenterPosition

◆ fConstructed

G4bool G4VScoringMesh::fConstructed
protected

Definition at line 225 of file G4VScoringMesh.hh.

Referenced by Construct(), and WorkerConstruct().

◆ fCurrentPS

◆ fDivisionAxisNames

G4String G4VScoringMesh::fDivisionAxisNames[3]
protected

◆ fDrawPSName

◆ fDrawUnit

◆ fDrawUnitValue

G4double G4VScoringMesh::fDrawUnitValue
protected

◆ fGeometryHasBeenDestroyed

G4bool G4VScoringMesh::fGeometryHasBeenDestroyed
protected

Definition at line 252 of file G4VScoringMesh.hh.

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

◆ fMap

◆ fMeshElementLogical

◆ fMFD

◆ fNSegment

◆ fParallelWorldProcess

G4ParallelWorldProcess* G4VScoringMesh::fParallelWorldProcess
protected

Definition at line 251 of file G4VScoringMesh.hh.

Referenced by GetParallelWorldProcess(), and SetParallelWorldProcess().

◆ fRotationMatrix

◆ fShape

◆ fSize

◆ fWorldName

◆ layeredMassFlg

G4bool G4VScoringMesh::layeredMassFlg
protected

Definition at line 259 of file G4VScoringMesh.hh.

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

◆ nMeshIsSet

G4bool G4VScoringMesh::nMeshIsSet
protected

Definition at line 241 of file G4VScoringMesh.hh.

Referenced by ReadyForQuantity(), and SetNumberOfSegments().

◆ sizeIsSet

G4bool G4VScoringMesh::sizeIsSet
protected

Definition at line 240 of file G4VScoringMesh.hh.

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

◆ verboseLevel


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