Geant4-11
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes
G4ScoringProbe Class Reference

#include <G4ScoringProbe.hh>

Inheritance diagram for G4ScoringProbe:
G4VScoringMesh

Public Types

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

Public Member Functions

void Accumulate (G4THitsMap< G4double > *map)
 
void Accumulate (G4THitsMap< G4StatDouble > *map)
 
void Activate (G4bool vl=true)
 
virtual void Construct (G4VPhysicalVolume *fWorldPhys)
 
virtual void Draw (RunScore *, G4VScoreColorMap *, G4int)
 
virtual void DrawColumn (RunScore *, G4VScoreColorMap *, G4int, G4int)
 
void DrawMesh (const G4String &psName, G4int idxPlane, G4int iColumn, G4VScoreColorMap *colorMap)
 
void DrawMesh (const G4String &psName, G4VScoreColorMap *colorMap, G4int axflg=111)
 
void Dump ()
 
G4bool FindPrimitiveScorer (const G4String &psname)
 
 G4ScoringProbe (G4String lvName, G4double half_size, G4bool checkOverlap=false)
 
void GeometryHasBeenDestroyed ()
 
G4double GetAngleSpan () const
 
G4int GetCopyNumberLevel () const
 
G4String GetCurrentPSUnit ()
 
void GetDivisionAxisNames (G4String divisionAxisNames[3])
 
G4LogicalVolumeGetMeshElementLogical () const
 
G4int GetNumberOfProbes () const
 
void GetNumberOfSegments (G4int nSegment[3])
 
G4ParallelWorldProcessGetParallelWorldProcess () const
 
G4VPrimitiveScorerGetPrimitiveScorer (const G4String &name)
 
G4double GetProbeSize () const
 
G4String GetPSUnit (const G4String &psname)
 
G4double GetPSUnitValue (const G4String &psname)
 
G4RotationMatrix GetRotationMatrix () const
 
MeshScoreMap GetScoreMap () const
 
MeshShape GetShape () const
 
G4ThreeVector GetSize () const
 
G4double GetStartAngle () const
 
G4ThreeVector GetTranslation () const
 
const G4StringGetWorldName () const
 
G4bool IsActive () const
 
G4bool IsCurrentPrimitiveScorerNull ()
 
G4bool LayeredMassFlg ()
 
virtual void List () const
 
void LocateProbe (G4ThreeVector pos)
 
void Merge (const G4VScoringMesh *scMesh)
 
G4bool ReadyForQuantity () const
 
void ResetScore ()
 
void RotateX (G4double delta)
 
void RotateY (G4double delta)
 
void RotateZ (G4double delta)
 
void SetAngles (G4double, G4double)
 
void SetCenterPosition (G4double centerPosition[3])
 
void SetCopyNumberLevel (G4int val)
 
void SetCurrentPrimitiveScorer (const G4String &name)
 
void SetCurrentPSUnit (const G4String &unit)
 
void SetDrawPSName (const G4String &psname)
 
void SetFilter (G4VSDFilter *filter)
 
G4bool SetMaterial (G4String val)
 
void SetMeshElementLogical (G4LogicalVolume *val)
 
void SetNullToCurrentPrimitiveScorer ()
 
void SetNumberOfSegments (G4int nSegment[3])
 
void SetParallelWorldProcess (G4ParallelWorldProcess *proc)
 
void SetPrimitiveScorer (G4VPrimitiveScorer *ps)
 
void SetProbeSize (G4double val)
 
void SetSize (G4double size[3])
 
void SetVerboseLevel (G4int vl)
 
virtual void WorkerConstruct (G4VPhysicalVolume *fWorldPhys)
 
 ~G4ScoringProbe ()
 

Protected Member Functions

virtual void SetupGeometry (G4VPhysicalVolume *)
 

Protected Attributes

G4bool chkOverlap
 
G4int copyNumberLevel
 
G4bool fActive
 
G4double fAngle [2]
 
G4ThreeVector fCenterPosition
 
G4bool fConstructed
 
G4VPrimitiveScorerfCurrentPS
 
G4String fDivisionAxisNames [3]
 
G4String fDrawPSName
 
G4String fDrawUnit
 
G4double fDrawUnitValue
 
G4bool fGeometryHasBeenDestroyed
 
MeshScoreMap fMap
 
G4LogicalVolumefMeshElementLogical
 
G4MultiFunctionalDetectorfMFD
 
G4int fNSegment [3]
 
G4ParallelWorldProcessfParallelWorldProcess
 
G4RotationMatrixfRotationMatrix
 
MeshShape fShape
 
G4double fSize [3]
 
G4String fWorldName
 
G4bool layeredMassFlg
 
G4MateriallayeredMaterial
 
G4String layeredMaterialName
 
G4String logVolName
 
G4bool nMeshIsSet
 
std::vector< G4ThreeVectorposVec
 
G4double probeSize
 
G4String regName
 
G4bool sizeIsSet
 
G4int verboseLevel
 

Detailed Description

Definition at line 37 of file G4ScoringProbe.hh.

Member Typedef Documentation

◆ EventScore

Definition at line 66 of file G4VScoringMesh.hh.

◆ MeshScoreMap

using G4VScoringMesh::MeshScoreMap = std::map<G4String, RunScore*>
inherited

Definition at line 68 of file G4VScoringMesh.hh.

◆ RunScore

Definition at line 67 of file G4VScoringMesh.hh.

Member Enumeration Documentation

◆ MeshShape

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

Definition at line 57 of file G4VScoringMesh.hh.

58 {
59 box,
60 cylinder,
61 sphere,
62 realWorldLogVol,
63 probe,
64 undefined = -1
65 };

Constructor & Destructor Documentation

◆ G4ScoringProbe()

G4ScoringProbe::G4ScoringProbe ( G4String  lvName,
G4double  half_size,
G4bool  checkOverlap = false 
)

Definition at line 50 of file G4ScoringProbe.cc.

52 : G4VScoringMesh(lvName)
53 , chkOverlap(checkOverlap)
54 , layeredMaterialName("none")
55 , layeredMaterial(nullptr)
56{
58 logVolName = lvName;
59 probeSize = half_size;
60 G4double hs[] = { half_size, half_size, half_size };
61 SetSize(hs);
62 G4int nBin[] = { 1, 1, 1 };
64 regName = lvName + "_region";
66 {
67 new G4Region(regName);
68 }
69}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4double probeSize
G4String logVolName
G4String layeredMaterialName
G4Material * layeredMaterial
void SetNumberOfSegments(G4int nSegment[3])
G4VScoringMesh(const G4String &wName)
void SetSize(G4double size[3])
G4bool IsMasterThread()
Definition: G4Threading.cc:124

References G4VScoringMesh::fShape, G4Threading::IsMasterThread(), logVolName, G4VScoringMesh::probe, probeSize, regName, G4VScoringMesh::SetNumberOfSegments(), and G4VScoringMesh::SetSize().

◆ ~G4ScoringProbe()

G4ScoringProbe::~G4ScoringProbe ( )

Definition at line 71 of file G4ScoringProbe.cc.

71{}

Member Function Documentation

◆ Accumulate() [1/2]

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

Definition at line 386 of file G4VScoringMesh.cc.

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

References G4VScoringMesh::fMap, G4cout, G4endl, anonymous_namespace{G4QuasiElRatios.cc}::map, and G4VScoringMesh::verboseLevel.

◆ Accumulate() [2/2]

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

Definition at line 410 of file G4VScoringMesh.cc.

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

References G4VScoringMesh::fMap, G4cout, G4endl, anonymous_namespace{G4QuasiElRatios.cc}::map, and G4VScoringMesh::verboseLevel.

◆ Activate()

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

Definition at line 95 of file G4VScoringMesh.hh.

95{ fActive = vl; }

References G4VScoringMesh::fActive.

◆ Construct()

void G4VScoringMesh::Construct ( G4VPhysicalVolume fWorldPhys)
virtualinherited

Definition at line 434 of file G4VScoringMesh.cc.

435{
436 if(fConstructed)
437 {
439 {
440 SetupGeometry(fWorldPhys);
442 }
443 if(verboseLevel > 0)
444 G4cout << fWorldName << " --- All quantities are reset." << G4endl;
445 ResetScore();
446 }
447 else
448 {
449 fConstructed = true;
450 SetupGeometry(fWorldPhys);
451 }
452}
virtual void SetupGeometry(G4VPhysicalVolume *fWorldPhys)=0
G4bool fGeometryHasBeenDestroyed

References G4VScoringMesh::fConstructed, G4VScoringMesh::fGeometryHasBeenDestroyed, G4VScoringMesh::fWorldName, G4cout, G4endl, G4VScoringMesh::ResetScore(), G4VScoringMesh::SetupGeometry(), and G4VScoringMesh::verboseLevel.

Referenced by G4RunManager::ConstructScoringWorlds().

◆ Draw()

virtual void G4ScoringProbe::Draw ( RunScore ,
G4VScoreColorMap ,
G4int   
)
inlinevirtual

Implements G4VScoringMesh.

Definition at line 74 of file G4ScoringProbe.hh.

76 {
77 ;
78 }

◆ DrawColumn()

virtual void G4ScoringProbe::DrawColumn ( RunScore ,
G4VScoreColorMap ,
G4int  ,
G4int   
)
inlinevirtual

Implements G4VScoringMesh.

Definition at line 79 of file G4ScoringProbe.hh.

81 {
82 ;
83 }

◆ DrawMesh() [1/2]

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

Definition at line 368 of file G4VScoringMesh.cc.

370{
371 fDrawPSName = psName;
372 MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
373 if(fMapItr != fMap.end())
374 {
375 fDrawUnit = GetPSUnit(psName);
377 DrawColumn(fMapItr->second, colorMap, idxPlane, iColumn);
378 }
379 else
380 {
381 G4cerr << "Scorer <" << psName << "> is not defined. Method ignored."
382 << G4endl;
383 }
384}
G4GLOB_DLL std::ostream G4cerr
G4double fDrawUnitValue
G4String fDrawPSName
G4double GetPSUnitValue(const G4String &psname)
G4String GetPSUnit(const G4String &psname)
virtual void DrawColumn(RunScore *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)=0

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

◆ DrawMesh() [2/2]

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

Definition at line 350 of file G4VScoringMesh.cc.

352{
353 fDrawPSName = psName;
354 MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
355 if(fMapItr != fMap.end())
356 {
357 fDrawUnit = GetPSUnit(psName);
359 Draw(fMapItr->second, colorMap, axflg);
360 }
361 else
362 {
363 G4cerr << "Scorer <" << psName << "> is not defined. Method ignored."
364 << G4endl;
365 }
366}
virtual void Draw(RunScore *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0

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

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

◆ Dump()

void G4VScoringMesh::Dump ( )
inherited

Definition at line 338 of file G4VScoringMesh.cc.

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

References G4VScoringMesh::fMap, G4VScoringMesh::fWorldName, G4cout, and G4endl.

◆ FindPrimitiveScorer()

G4bool G4VScoringMesh::FindPrimitiveScorer ( const G4String psname)
inherited

Definition at line 223 of file G4VScoringMesh.cc.

224{
225 MeshScoreMap::iterator itr = fMap.find(psname);
226 if(itr == fMap.end())
227 return false;
228 return true;
229}

References G4VScoringMesh::fMap.

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

◆ GeometryHasBeenDestroyed()

void G4VScoringMesh::GeometryHasBeenDestroyed ( )
inlineinherited

◆ GetAngleSpan()

G4double G4VScoringMesh::GetAngleSpan ( ) const
inlineinherited

Definition at line 129 of file G4VScoringMesh.hh.

129{ return fAngle[1]; }
G4double fAngle[2]

References G4VScoringMesh::fAngle.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetCopyNumberLevel()

G4int G4VScoringMesh::GetCopyNumberLevel ( ) const
inlineinherited

Definition at line 263 of file G4VScoringMesh.hh.

263{ return copyNumberLevel; }

References G4VScoringMesh::copyNumberLevel.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetCurrentPSUnit()

G4String G4VScoringMesh::GetCurrentPSUnit ( )
inherited

Definition at line 244 of file G4VScoringMesh.cc.

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

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

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetDivisionAxisNames()

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

Definition at line 287 of file G4VScoringMesh.cc.

288{
289 for(int i = 0; i < 3; i++)
290 divisionAxisNames[i] = fDivisionAxisNames[i];
291}
G4String fDivisionAxisNames[3]

References G4VScoringMesh::fDivisionAxisNames.

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

◆ GetMeshElementLogical()

G4LogicalVolume * G4VScoringMesh::GetMeshElementLogical ( ) const
inlineinherited

Definition at line 232 of file G4VScoringMesh.hh.

233 {
234 return fMeshElementLogical;
235 }

References G4VScoringMesh::fMeshElementLogical.

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

◆ GetNumberOfProbes()

G4int G4ScoringProbe::GetNumberOfProbes ( ) const
inline

Definition at line 64 of file G4ScoringProbe.hh.

64{ return posVec.size(); }
std::vector< G4ThreeVector > posVec

References posVec.

◆ GetNumberOfSegments()

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

Definition at line 140 of file G4VScoringMesh.cc.

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

References G4VScoringMesh::fNSegment.

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

◆ GetParallelWorldProcess()

G4ParallelWorldProcess * G4VScoringMesh::GetParallelWorldProcess ( ) const
inlineinherited

Definition at line 246 of file G4VScoringMesh.hh.

247 {
249 }
G4ParallelWorldProcess * fParallelWorldProcess

References G4VScoringMesh::fParallelWorldProcess.

Referenced by G4RunManager::ConstructScoringWorlds(), and G4WorkerRunManager::ConstructScoringWorlds().

◆ GetPrimitiveScorer()

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

◆ GetProbeSize()

G4double G4ScoringProbe::GetProbeSize ( ) const
inline

Definition at line 66 of file G4ScoringProbe.hh.

66{ return probeSize; }

References probeSize.

◆ GetPSUnit()

G4String G4VScoringMesh::GetPSUnit ( const G4String psname)
inherited

Definition at line 231 of file G4VScoringMesh.cc.

232{
233 MeshScoreMap::iterator itr = fMap.find(psname);
234 if(itr == fMap.end())
235 {
236 return G4String("");
237 }
238 else
239 {
240 return GetPrimitiveScorer(psname)->GetUnit();
241 }
242}
G4VPrimitiveScorer * GetPrimitiveScorer(const G4String &name)

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

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

◆ GetPSUnitValue()

G4double G4VScoringMesh::GetPSUnitValue ( const G4String psname)
inherited

Definition at line 274 of file G4VScoringMesh.cc.

275{
276 MeshScoreMap::iterator itr = fMap.find(psname);
277 if(itr == fMap.end())
278 {
279 return 1.;
280 }
281 else
282 {
283 return GetPrimitiveScorer(psname)->GetUnitValue();
284 }
285}
G4double GetUnitValue() const

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

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

◆ GetRotationMatrix()

G4RotationMatrix G4VScoringMesh::GetRotationMatrix ( ) const
inlineinherited

Definition at line 141 of file G4VScoringMesh.hh.

142 {
144 return *fRotationMatrix;
145 else
147 }
static DLL_API const HepRotation IDENTITY
Definition: Rotation.h:366
G4RotationMatrix * fRotationMatrix

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

Referenced by G4GMocrenFileSceneHandler::AddSolid().

◆ GetScoreMap()

MeshScoreMap G4VScoringMesh::GetScoreMap ( ) const
inlineinherited

◆ GetShape()

MeshShape G4VScoringMesh::GetShape ( ) const
inlineinherited

◆ GetSize()

G4ThreeVector G4VScoringMesh::GetSize ( ) const
inherited

Definition at line 105 of file G4VScoringMesh.cc.

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

References G4VScoringMesh::fSize, and G4VScoringMesh::sizeIsSet.

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

◆ GetStartAngle()

G4double G4VScoringMesh::GetStartAngle ( ) const
inlineinherited

Definition at line 128 of file G4VScoringMesh.hh.

128{ return fAngle[0]; }

References G4VScoringMesh::fAngle.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetTranslation()

G4ThreeVector G4VScoringMesh::GetTranslation ( ) const
inlineinherited

Definition at line 133 of file G4VScoringMesh.hh.

133{ return fCenterPosition; }
G4ThreeVector fCenterPosition

References G4VScoringMesh::fCenterPosition.

Referenced by G4GMocrenFileSceneHandler::AddSolid().

◆ GetWorldName()

const G4String & G4VScoringMesh::GetWorldName ( ) const
inlineinherited

◆ IsActive()

G4bool G4VScoringMesh::IsActive ( ) const
inlineinherited

Definition at line 93 of file G4VScoringMesh.hh.

93{ return fActive; }

References G4VScoringMesh::fActive.

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

◆ IsCurrentPrimitiveScorerNull()

G4bool G4VScoringMesh::IsCurrentPrimitiveScorerNull ( )
inlineinherited

Definition at line 164 of file G4VScoringMesh.hh.

165 {
166 if(fCurrentPS == nullptr)
167 return true;
168 else
169 return false;
170 }

References G4VScoringMesh::fCurrentPS.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ LayeredMassFlg()

G4bool G4VScoringMesh::LayeredMassFlg ( )
inlineinherited

◆ List()

void G4ScoringProbe::List ( ) const
virtual

Reimplemented from G4VScoringMesh.

Definition at line 73 of file G4ScoringProbe.cc.

74{
75 G4cout << "G4ScoringProbe : " << logVolName << G4endl;
76 G4int np = posVec.size();
77 for(G4int i = 0; i < np; i++)
78 {
79 G4cout << " >> probe #" << i << " at " << posVec[i] << G4endl;
80 }
82}
virtual void List() const

References G4cout, G4endl, G4VScoringMesh::List(), logVolName, and posVec.

◆ LocateProbe()

void G4ScoringProbe::LocateProbe ( G4ThreeVector  pos)
inline

Definition at line 58 of file G4ScoringProbe.hh.

59 {
60 posVec.push_back(pos);
61 G4int nbin[] = { static_cast<G4int>(posVec.size()), 1, 1 };
63 }
static const G4double pos

References pos, posVec, and G4VScoringMesh::SetNumberOfSegments().

◆ Merge()

void G4VScoringMesh::Merge ( const G4VScoringMesh scMesh)
inherited

Definition at line 476 of file G4VScoringMesh.cc.

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

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

Referenced by G4ScoringManager::Merge().

◆ ReadyForQuantity()

G4bool G4VScoringMesh::ReadyForQuantity ( ) const
inlineinherited

Definition at line 192 of file G4VScoringMesh.hh.

192{ return (sizeIsSet && nMeshIsSet); }

References G4VScoringMesh::nMeshIsSet, and G4VScoringMesh::sizeIsSet.

Referenced by G4VScoringMesh::SetPrimitiveScorer().

◆ ResetScore()

void G4VScoringMesh::ResetScore ( )
inherited

Definition at line 75 of file G4VScoringMesh.cc.

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

References G4VScoringMesh::fMap, G4cout, G4endl, and G4VScoringMesh::verboseLevel.

Referenced by G4VScoringMesh::Construct(), and G4VScoringMesh::WorkerConstruct().

◆ RotateX()

void G4VScoringMesh::RotateX ( G4double  delta)
inherited

Definition at line 145 of file G4VScoringMesh.cc.

146{
147 if(!fRotationMatrix)
149 fRotationMatrix->rotateX(delta);
150}
CLHEP::HepRotation G4RotationMatrix
HepRotation & rotateX(double delta)
Definition: Rotation.cc:61

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

Referenced by G4ScoringMessenger::SetNewValue().

◆ RotateY()

void G4VScoringMesh::RotateY ( G4double  delta)
inherited

Definition at line 152 of file G4VScoringMesh.cc.

153{
154 if(!fRotationMatrix)
156 fRotationMatrix->rotateY(delta);
157}
HepRotation & rotateY(double delta)
Definition: Rotation.cc:74

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

Referenced by G4ScoringMessenger::SetNewValue().

◆ RotateZ()

void G4VScoringMesh::RotateZ ( G4double  delta)
inherited

Definition at line 159 of file G4VScoringMesh.cc.

160{
161 if(!fRotationMatrix)
163 fRotationMatrix->rotateZ(delta);
164}
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:87

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

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetAngles()

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

Definition at line 112 of file G4VScoringMesh.cc.

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

References G4VScoringMesh::fAngle.

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetCenterPosition()

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

Definition at line 118 of file G4VScoringMesh.cc.

119{
121 G4ThreeVector(centerPosition[0], centerPosition[1], centerPosition[2]);
122}

References G4VScoringMesh::fCenterPosition.

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetCopyNumberLevel()

void G4VScoringMesh::SetCopyNumberLevel ( G4int  val)
inlineinherited

Definition at line 262 of file G4VScoringMesh.hh.

262{ copyNumberLevel = val; }

References G4VScoringMesh::copyNumberLevel.

◆ SetCurrentPrimitiveScorer()

void G4VScoringMesh::SetCurrentPrimitiveScorer ( const G4String name)
inherited

Definition at line 212 of file G4VScoringMesh.cc.

213{
215 if(!fCurrentPS)
216 {
217 G4cerr << "ERROR : G4VScoringMesh::SetCurrentPrimitiveScorer() : The "
218 "primitive scorer <"
219 << name << "> does not found." << G4endl;
220 }
221}

References G4VScoringMesh::fCurrentPS, G4cerr, G4endl, G4VScoringMesh::GetPrimitiveScorer(), and G4InuclParticleNames::name().

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetCurrentPSUnit()

void G4VScoringMesh::SetCurrentPSUnit ( const G4String unit)
inherited

Definition at line 260 of file G4VScoringMesh.cc.

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

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

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetDrawPSName()

void G4VScoringMesh::SetDrawPSName ( const G4String psname)
inlineinherited

Definition at line 180 of file G4VScoringMesh.hh.

180{ fDrawPSName = psname; }

References G4VScoringMesh::fDrawPSName.

◆ SetFilter()

void G4VScoringMesh::SetFilter ( G4VSDFilter filter)
inherited

Definition at line 190 of file G4VScoringMesh.cc.

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

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

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

◆ SetMaterial()

G4bool G4ScoringProbe::SetMaterial ( G4String  val)

Definition at line 131 of file G4ScoringProbe.cc.

132{
134 {
135 if(val == "none")
136 {
138 layeredMassFlg = false;
139 layeredMaterial = nullptr;
140 }
141 else
142 {
144 if(!mat)
145 {
146 return false;
147 }
149 layeredMassFlg = true;
150 layeredMaterial = mat;
151 }
153 assert(region != nullptr);
154 region->UpdateMaterialList();
155 }
156 return true;
157}
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const

References G4NistManager::FindOrBuildMaterial(), G4RegionStore::GetInstance(), G4RegionStore::GetRegion(), G4NistManager::Instance(), G4Threading::IsMasterThread(), G4VScoringMesh::layeredMassFlg, layeredMaterial, layeredMaterialName, and regName.

◆ SetMeshElementLogical()

void G4VScoringMesh::SetMeshElementLogical ( G4LogicalVolume val)
inlineinherited

Definition at line 228 of file G4VScoringMesh.hh.

229 {
231 }

References G4VScoringMesh::fMeshElementLogical.

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

◆ SetNullToCurrentPrimitiveScorer()

void G4VScoringMesh::SetNullToCurrentPrimitiveScorer ( )
inlineinherited

Definition at line 186 of file G4VScoringMesh.hh.

186{ fCurrentPS = nullptr; }

References G4VScoringMesh::fCurrentPS.

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

◆ SetNumberOfSegments()

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

Definition at line 123 of file G4VScoringMesh.cc.

124{
127 {
128 for(int 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)
Definition: G4Exception.cc:35

References G4VScoringMesh::fNSegment, G4VScoringMesh::fShape, G4Exception(), JustWarning, G4VScoringMesh::nMeshIsSet, G4VScoringMesh::probe, and G4VScoringMesh::realWorldLogVol.

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

◆ SetParallelWorldProcess()

void G4VScoringMesh::SetParallelWorldProcess ( G4ParallelWorldProcess proc)
inlineinherited

◆ SetPrimitiveScorer()

void G4VScoringMesh::SetPrimitiveScorer ( G4VPrimitiveScorer ps)
inherited

Definition at line 166 of file G4VScoringMesh.cc.

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

References G4VScoringMesh::fCurrentPS, G4VScoringMesh::fMap, G4VScoringMesh::fMFD, G4VScoringMesh::fNSegment, G4VScoringMesh::fWorldName, G4cerr, G4cout, G4endl, G4VPrimitiveScorer::GetName(), anonymous_namespace{G4QuasiElRatios.cc}::map, G4VScoringMesh::ReadyForQuantity(), G4MultiFunctionalDetector::RegisterPrimitive(), G4VPrimitiveScorer::SetNijk(), and G4VScoringMesh::verboseLevel.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetProbeSize()

void G4ScoringProbe::SetProbeSize ( G4double  val)
inline

Definition at line 65 of file G4ScoringProbe.hh.

65{ probeSize = val; }

References probeSize.

◆ SetSize()

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

Definition at line 87 of file G4VScoringMesh.cc.

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

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

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

◆ SetupGeometry()

void G4ScoringProbe::SetupGeometry ( G4VPhysicalVolume worldPhys)
protectedvirtual

Implements G4VScoringMesh.

Definition at line 90 of file G4ScoringProbe.cc.

91{
93 {
94 auto worldLog = worldPhys->GetLogicalVolume();
96 assert(region != nullptr);
97 region->AddRootLogicalVolume(worldLog);
98 region->SetWorld(worldPhys);
99
100 auto boxSolid =
101 new G4Box(logVolName + "_solid", probeSize, probeSize, probeSize);
103 new G4LogicalVolume(boxSolid, layeredMaterial, logVolName + "_log");
104
105 G4int np = posVec.size();
106 for(G4int i = 0; i < np; i++)
107 {
109 worldLog, false, i, chkOverlap);
110 }
111
112 G4VisAttributes* wisatt = new G4VisAttributes(G4Colour(.5, .5, .5));
113 wisatt->SetVisibility(false);
114 worldLog->SetVisAttributes(wisatt);
115 G4VisAttributes* visatt = new G4VisAttributes(G4Colour(.5, .5, .5));
116 visatt->SetVisibility(true);
118 }
119 else
120 {
124 assert(fMeshElementLogical != nullptr);
125 l.unlock();
126 }
127
129}
Definition: G4Box.hh:56
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true, G4bool reverseSearch=false) const
static G4LogicalVolumeStore * GetInstance()
void SetVisAttributes(const G4VisAttributes *pVA)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
G4LogicalVolume * GetLogicalVolume() const
void SetVisibility(G4bool=true)

References chkOverlap, G4VScoringMesh::fMeshElementLogical, G4VScoringMesh::fMFD, G4LogicalVolumeStore::GetInstance(), G4RegionStore::GetInstance(), G4VPhysicalVolume::GetLogicalVolume(), G4RegionStore::GetRegion(), G4LogicalVolumeStore::GetVolume(), G4Threading::IsMasterThread(), layeredMaterial, anonymous_namespace{G4ScoringProbe.cc}::logvolmutex, logVolName, posVec, probeSize, regName, G4LogicalVolume::SetSensitiveDetector(), G4LogicalVolume::SetVisAttributes(), G4VisAttributes::SetVisibility(), and G4TemplateAutoLock< _Mutex_t >::unlock().

◆ SetVerboseLevel()

void G4VScoringMesh::SetVerboseLevel ( G4int  vl)
inlineinherited

Definition at line 188 of file G4VScoringMesh.hh.

188{ verboseLevel = vl; }

References G4VScoringMesh::verboseLevel.

Referenced by G4ScoringManager::RegisterScoringMesh().

◆ WorkerConstruct()

void G4VScoringMesh::WorkerConstruct ( G4VPhysicalVolume fWorldPhys)
virtualinherited

Field Documentation

◆ chkOverlap

G4bool G4ScoringProbe::chkOverlap
protected

Definition at line 52 of file G4ScoringProbe.hh.

Referenced by SetupGeometry().

◆ copyNumberLevel

G4int G4VScoringMesh::copyNumberLevel
protectedinherited

◆ fActive

G4bool G4VScoringMesh::fActive
protectedinherited

Definition at line 202 of file G4VScoringMesh.hh.

Referenced by G4VScoringMesh::Activate(), and G4VScoringMesh::IsActive().

◆ fAngle

G4double G4VScoringMesh::fAngle[2]
protectedinherited

◆ fCenterPosition

G4ThreeVector G4VScoringMesh::fCenterPosition
protectedinherited

◆ fConstructed

G4bool G4VScoringMesh::fConstructed
protectedinherited

◆ fCurrentPS

G4VPrimitiveScorer* G4VScoringMesh::fCurrentPS
protectedinherited

◆ fDivisionAxisNames

G4String G4VScoringMesh::fDivisionAxisNames[3]
protectedinherited

◆ fDrawPSName

G4String G4VScoringMesh::fDrawPSName
protectedinherited

◆ fDrawUnit

G4String G4VScoringMesh::fDrawUnit
protectedinherited

◆ fDrawUnitValue

G4double G4VScoringMesh::fDrawUnitValue
protectedinherited

◆ fGeometryHasBeenDestroyed

G4bool G4VScoringMesh::fGeometryHasBeenDestroyed
protectedinherited

◆ fMap

MeshScoreMap G4VScoringMesh::fMap
protectedinherited

◆ fMeshElementLogical

G4LogicalVolume* G4VScoringMesh::fMeshElementLogical
protectedinherited

◆ fMFD

G4MultiFunctionalDetector* G4VScoringMesh::fMFD
protectedinherited

◆ fNSegment

G4int G4VScoringMesh::fNSegment[3]
protectedinherited

◆ fParallelWorldProcess

G4ParallelWorldProcess* G4VScoringMesh::fParallelWorldProcess
protectedinherited

◆ fRotationMatrix

G4RotationMatrix* G4VScoringMesh::fRotationMatrix
protectedinherited

◆ fShape

MeshShape G4VScoringMesh::fShape
protectedinherited

◆ fSize

G4double G4VScoringMesh::fSize[3]
protectedinherited

◆ fWorldName

G4String G4VScoringMesh::fWorldName
protectedinherited

◆ layeredMassFlg

G4bool G4VScoringMesh::layeredMassFlg
protectedinherited

Definition at line 269 of file G4VScoringMesh.hh.

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

◆ layeredMaterial

G4Material* G4ScoringProbe::layeredMaterial
protected

Definition at line 54 of file G4ScoringProbe.hh.

Referenced by SetMaterial(), and SetupGeometry().

◆ layeredMaterialName

G4String G4ScoringProbe::layeredMaterialName
protected

Definition at line 53 of file G4ScoringProbe.hh.

Referenced by SetMaterial().

◆ logVolName

G4String G4ScoringProbe::logVolName
protected

Definition at line 49 of file G4ScoringProbe.hh.

Referenced by G4ScoringProbe(), List(), and SetupGeometry().

◆ nMeshIsSet

G4bool G4VScoringMesh::nMeshIsSet
protectedinherited

◆ posVec

std::vector<G4ThreeVector> G4ScoringProbe::posVec
protected

Definition at line 50 of file G4ScoringProbe.hh.

Referenced by GetNumberOfProbes(), List(), LocateProbe(), and SetupGeometry().

◆ probeSize

G4double G4ScoringProbe::probeSize
protected

Definition at line 51 of file G4ScoringProbe.hh.

Referenced by G4ScoringProbe(), GetProbeSize(), SetProbeSize(), and SetupGeometry().

◆ regName

G4String G4ScoringProbe::regName
protected

Definition at line 55 of file G4ScoringProbe.hh.

Referenced by G4ScoringProbe(), SetMaterial(), and SetupGeometry().

◆ sizeIsSet

G4bool G4VScoringMesh::sizeIsSet
protectedinherited

◆ verboseLevel

G4int G4VScoringMesh::verboseLevel
protectedinherited

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