Hits¶
Hit¶
A hit is a snapshot of the physical interaction of a track in the
sensitive region of a detector. In it you can store information
associated with a G4Step
object. This information can be
the position and time of the step,
the momentum and energy of the track,
the energy deposition of the step,
geometrical information,
or any combination of the above.
G4VHit
¶
G4VHit
is an abstract base class which represents a hit. You must
inherit this base class and derive your own concrete hit class(es). The
member data of your concrete hit class can be, and should be, your
choice.
As with G4THitsCollection
, authors of subclasses must declare
templated G4Allocator
s for their class. They must also implement
operator new() and operator delete() which use these allocators.
G4VHit
has two virtual methods, Draw()
and Print()
. To draw
or print out your concrete hits, these methods should be implemented.
How to define the drawing method is described in
Polylines, Markers and Text.
G4THitsCollection
¶
G4VHit
is an abstract class from which you derive your own concrete
classes. During the processing of a given event, represented by a
G4Event
object, many objects of the hit class will be produced,
collected and associated with the event. Therefore, for each concrete
hit class you must also prepare a concrete class derived from
G4VHitsCollection
, an abstract class which represents a vector
collection of user defined hits.
G4THitsCollection
is a template class derived from
G4VHitsCollection
, and the concrete hit collection class of a
particular G4VHit
concrete class can be instantiated from this
template class. Each object of a hit collection must have a unique name
for each event.
G4Event
has a G4HCofThisEvent
class object, that is a container
class of collections of hits. Hit collections are stored by their
pointers, whose type is that of the base class.
An example of a concrete hit class¶
Listing 53 shows an example of a concrete hit class.
//============ header file =====================
#ifndef B2TrackerHit_h
#define B2TrackerHit_h 1
#include "G4VHit.hh"
#include "G4THitsCollection.hh"
#include "G4Allocator.hh"
#include "G4ThreeVector.hh"
#include "tls.hh"
namespace B2
{
/// Tracker hit class
///
/// It defines data members to store the trackID, chamberNb, energy deposit,
/// and position of charged particles in a selected volume:
/// - fTrackID, fChamberNB, fEdep, fPos
class TrackerHit : public G4VHit
{
public:
TrackerHit();
TrackerHit(const TrackerHit&) = default;
~TrackerHit() override;
// operators
TrackerHit& operator=(const TrackerHit&) = default;
G4bool operator==(const TrackerHit&) const;
inline void* operator new(size_t);
inline void operator delete(void*);
// methods from base class
void Draw() override;
void Print() override;
// Set methods
void SetTrackID (G4int track) { fTrackID = track; };
void SetChamberNb(G4int chamb) { fChamberNb = chamb; };
void SetEdep (G4double de) { fEdep = de; };
void SetPos (G4ThreeVector xyz){ fPos = xyz; };
// Get methods
G4int GetTrackID() const { return fTrackID; };
G4int GetChamberNb() const { return fChamberNb; };
G4double GetEdep() const { return fEdep; };
G4ThreeVector GetPos() const { return fPos; };
private:
G4int fTrackID = -1;
G4int fChamberNb = -1;
G4double fEdep = 0.;
G4ThreeVector fPos;
};
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
typedef G4THitsCollection<TrackerHit> TrackerHitsCollection;
extern G4ThreadLocal G4Allocator<TrackerHit>* TrackerHitAllocator;
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
inline void* TrackerHit::operator new(size_t)
{
if(!TrackerHitAllocator)
TrackerHitAllocator = new G4Allocator<TrackerHit>;
return (void *) TrackerHitAllocator->MallocSingle();
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
inline void TrackerHit::operator delete(void *hit)
{
TrackerHitAllocator->FreeSingle((TrackerHit*) hit);
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
}
#endif
//============ source file =====================
#include "TrackerHit.hh"
namespace B2
{
G4ThreadLocal G4Allocator<TrackerHit>* TrackerHit::TrackerHitAllocator = 0;
... snipped ...
G4Allocator
is a class for fast allocation of objects to the heap
through the paging mechanism. For details of G4Allocator
, refer to
General management classes. Use of G4Allocator
is not mandatory, but
it is recommended, especially for users who are not familiar with the
C++ memory allocation mechanism or alternative tools of memory
allocation. On the other hand, note that G4Allocator
is to be used
only for the concrete class that is not used as a base class of
any other classes. For example, do not use the G4Trajectory
class as a base class for a customized trajectory class, since
G4Trajectory
uses G4Allocator
.
G4THitsMap
¶
G4THitsMap
is an alternative to G4THitsCollection
.
G4THitsMap
does not demand G4VHit
, but instead any variable
which can be mapped with an integer key. Typically the key is a copy
number of the volume, and the mapped value could for example be a
double, such as the energy deposition in a volume. G4THitsMap
is
convenient for applications which do not need to output event-by-event
data but instead just accumulate them. All the G4VPrimitiveScorer
classes discussed in G4MultiFunctionalDetector and G4VPrimitiveScorer use G4THitsMap
.
G4THitsMap
is derived from the G4VHitsCollection
abstract base
class and all objects of this class are also stored in
G4HCofThisEvent
at the end of an event. How to access a
G4THitsMap
object is discussed in the following section
(G4MultiFunctionalDetector and G4VPrimitiveScorer).
Sensitive detector¶
G4VSensitiveDetector
¶
G4VSensitiveDetector
is an abstract base class which represents a
detector. The principal mandate of a sensitive detector is the
construction of hit objects using information from steps along a
particle track. The ProcessHits()
method of G4VSensitiveDetector
performs this task using G4Step
objects as input. The second
argument of ProcessHits()
method, i.e. G4TouchableHistory
, is
obsolete and not used. If user needs to define an artificial second
geometry, use Parallel Geometries.
ProcessHits()
method has a return type of G4bool
. This return
value is not used by Geant4 kernel. This return value may be used by
the user’s use-case where one sensitive detector dispatches ProcessHits()
to some subsequent (i.e. child) sensitive detectors, and to avoid
double-counting, one of these child detector may return true or false.
Your concrete detector class should be instantiated with the unique name of your detector. The name can be associated with one or more global names with “/” as a delimiter for categorizing your detectors. For example
myEMcal = new MyEMcal("/myDet/myCal/myEMcal");
where myEMcal
is the name of your detector. The detector must be
constructed in G4VUserDetectorConstruction::ConstructSDandField()
method. It must be assigned to one or more G4LogicalVolume
objects
to set the sensitivity of these volumes. Such assignment must be made in
the same G4VUserDetectorConstruction::ConstructSDandField()
method.
The pointer should also be registered to G4SDManager
, as described
in G4SDManager.
G4VSensitiveDetector
has three major virtual methods.
ProcessHits()
This method is invoked by
G4SteppingManager
when a step is composed in theG4LogicalVolume
which has the pointer to this sensitive detector. The first argument of this method is aG4Step
object of the current step. The second argument is aG4TouchableHistory
object for theReadout geometry
described in the next section. The second argument isNULL
ifReadout geometry
is not assigned to this sensitive detector. In this method, one or moreG4VHit
objects should be constructed if the current step is meaningful for your detector.Initialize()
This method is invoked at the beginning of each event. The argument of this method is an object of the
G4HCofThisEvent
class. Hit collections, where hits produced in this particular event are stored, can be associated with theG4HCofThisEvent
object in this method. The hit collections associated with theG4HCofThisEvent
object during this method can be used forduring the event processing
digitization.EndOfEvent()
This method is invoked at the end of each event. The argument of this method is the same object as the previous method. Hit collections occasionally created in your sensitive detector can be associated with the
G4HCofThisEvent
object.
G4SDManager¶
G4SDManager
is the singleton manager class for sensitive detectors.
Activation/inactivation of sensitive detectors
¶
The user interface commands activate
and inactivate
are
available to control your sensitive detectors. For example:
/hits/activate detector_name
/hits/inactivate detector_name
where detector_name
can be the detector name or the category name.
For example, if your EM calorimeter is named
/myDet/myCal/myEMcal
/hits/inactivate myCal
will inactivate all detectors belonging to the myCal
category.
Access to the hit collections
¶
Hit collections are accessed for various cases.
Digitization
Event filtering in
G4VUserStackingAction
“End of event” simple analysis
Drawing / printing hits
The following is an example of how to access the hit collection of a particular concrete type:
G4SDManager* fSDM = G4SDManager::GetSDMpointer();
G4RunManager* fRM = G4RunManager::GetRunManager();
G4int collectionID = fSDM->GetCollectionID("collection_name");
const G4Event* currentEvent = fRM->GetCurrentEvent();
G4HCofThisEvent* HCofEvent = currentEvent->GetHCofThisEvent();
MyHitsCollection* myCollection = (MyHitsCollection*)(HC0fEvent->GetHC(collectionID));
G4MultiFunctionalDetector
and G4VPrimitiveScorer
¶
G4MultiFunctionalDetector
is a concrete class derived from
G4VSensitiveDetector
. Instead of implementing a user-specific
detector class, G4MultiFunctionalDetector
allows the user to
register G4VPrimitiveScorer
classes to build up the sensitivity.
G4MultiFunctionalDetector
should be instantiated in the users
detector construction with its unique name and should be assigned to one
or more G4LogicalVolume
s.
G4VPrimitiveScorer
is an abstract base class representing a class to
be registered to G4MultiFunctionalDetector
that creates a
G4THitsMap
object of one physics quantity for an event. Geant4
provides many concrete primitive scorer classes listed in
Concrete classes of G4VPrimitiveScorer, and the user can also implement
his/her own primitive scorers. Each primitive scorer object must be
instantiated with a name that must be unique among primitive scorers
registered in a G4MultiFunctionalDetector
. Please note that a
primitive scorer object must not be shared by more than one
G4MultiFunctionalDetector
object.
As mentioned in Hit, each G4VPrimitiveScorer
generates one G4THitsMap
object per event. The name of the map
object is the same as the name of the primitive scorer. Each of the
concrete primitive scorers listed in
Concrete classes of G4VPrimitiveScorer generates a G4THitsMap<G4double>
that maps
a G4double
value to its key integer number. By default, the key is
taken as the copy number of the G4LogicalVolume
to which
G4MultiFunctionalDetector
is assigned. In case the logical volume is
uniquely placed in its mother volume and the mother is replicated, the
copy number of its mother volume can be taken by setting the second
argument of the G4VPrimitiveScorer
constructor, “depth” to 1, i.e.
one level up. Furthermore, in case the key must consider more than one
copy number of a different geometry hierarchy, the user can derive
his/her own primitive scorer from the provided concrete class and
implement the GetIndex(G4Step*)
virtual method to return the unique
key.
Listing 54 shows an example of primitive sensitivity class definitions.
void RE06DetectorConstruction::SetupDetectors()
{
G4String filterName, particleName;
G4SDParticleFilter* gammaFilter =
new G4SDParticleFilter(filterName="gammaFilter",particleName="gamma");
G4SDParticleFilter* electronFilter =
new G4SDParticleFilter(filterName="electronFilter",particleName="e-");
G4SDParticleFilter* positronFilter =
new G4SDParticleFilter(filterName="positronFilter",particleName="e+");
G4SDParticleFilter* epFilter = new G4SDParticleFilter(filterName="epFilter");
epFilter->add(particleName="e-");
epFilter->add(particleName="e+");
for(G4int i=0;i<3;i++)
{
for(G4int j=0;j<2;j++)
{
// Loop counter j = 0 : absorber
// = 1 : gap
G4String detName = fCalName[i];
if(j==0)
{ detName += "_abs"; }
else
{ detName += "_gap"; }
G4MultiFunctionalDetector* det = new G4MultiFunctionalDetector(detName);
// The second argument in each primitive means the "level" of geometrical hierarchy,
// the copy number of that level is used as the key of the G4THitsMap.
// For absorber (j = 0), the copy number of its own physical volume is used.
// For gap (j = 1), the copy number of its mother physical volume is used, since there
// is only one physical volume of gap is placed with respect to its mother.
G4VPrimitiveScorer* primitive;
primitive = new G4PSEnergyDeposit("eDep",j);
det->RegisterPrimitive(primitive);
primitive = new G4PSNofSecondary("nGamma",j);
primitive->SetFilter(gammaFilter);
det->RegisterPrimitive(primitive);
primitive = new G4PSNofSecondary("nElectron",j);
primitive->SetFilter(electronFilter);
det->RegisterPrimitive(primitive);
primitive = new G4PSNofSecondary("nPositron",j);
primitive->SetFilter(positronFilter);
det->RegisterPrimitive(primitive);
primitive = new G4PSMinKinEAtGeneration("minEkinGamma",j);
primitive->SetFilter(gammaFilter);
det->RegisterPrimitive(primitive);
primitive = new G4PSMinKinEAtGeneration("minEkinElectron",j);
primitive->SetFilter(electronFilter);
det->RegisterPrimitive(primitive);
primitive = new G4PSMinKinEAtGeneration("minEkinPositron",j);
primitive->SetFilter(positronFilter);
det->RegisterPrimitive(primitive);
primitive = new G4PSTrackLength("trackLength",j);
primitive->SetFilter(epFilter);
det->RegisterPrimitive(primitive);
primitive = new G4PSNofStep("nStep",j);
primitive->SetFilter(epFilter);
det->RegisterPrimitive(primitive);
G4SDManager::GetSDMpointer()->AddNewDetector(det);
if(j==0)
{ layerLogical[i]->SetSensitiveDetector(det); }
else
{ gapLogical[i]->SetSensitiveDetector(det); }
}
}
}
Each G4THitsMap
object can be accessed from G4HCofThisEvent
with
a unique collection ID number. This ID number can be obtained from
G4SDManager::GetCollectionID()
with a name of
G4MultiFunctionalDetector
and G4VPrimitiveScorer
connected with
a slush (“/”). G4THitsMap
has a [] operator taking the key value as
an argument and returning the pointer of the value. Please note that
the [] operator returns the pointer of the value. If you get zero
from the [] operator, it does not mean the value is zero, but that
the provided key does not exist. The value itself is accessible with an
asterisk (“*”). It is advised to check the validity of the returned
pointer before accessing the value. G4THitsMap
also has a +=
operator in order to accumulate event data into run data.
Listing 55 shows the use of
G4THitsMap
.
#include "Run.hh"
#include "G4RunManager.hh"
#include "G4Event.hh"
#include "G4SDManager.hh"
#include "G4HCofThisEvent.hh"
#include "G4THitsMap.hh"
#include "G4SystemOfUnits.hh"
namespace B3b
{
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Run::Run()
{ }
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Run::~Run()
{ }
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void Run::RecordEvent(const G4Event* event)
{
if ( fCollID_cryst < 0 ) {
fCollID_cryst
= G4SDManager::GetSDMpointer()->GetCollectionID("crystal/edep");
//G4cout << " fCollID_cryst: " << fCollID_cryst << G4endl;
}
if ( fCollID_patient < 0 ) {
fCollID_patient
= G4SDManager::GetSDMpointer()->GetCollectionID("patient/dose");
//G4cout << " fCollID_patient: " << fCollID_patient << G4endl;
}
G4int evtNb = event->GetEventID();
if (evtNb%fPrintModulo == 0) {
G4cout << G4endl << "---> end of event: " << evtNb << G4endl;
}
//Hits collections
//
G4HCofThisEvent* HCE = event->GetHCofThisEvent();
if(!HCE) return;
//Energy in crystals : identify 'good events'
//
const G4double eThreshold = 500*keV;
G4int nbOfFired = 0;
G4THitsMap<G4double>* evtMap =
static_cast<G4THitsMap<G4double>*>(HCE->GetHC(fCollID_cryst));
std::map<G4int,G4double*>::iterator itr;
for (itr = evtMap->GetMap()->begin(); itr != evtMap->GetMap()->end(); itr++) {
G4double edep = *(itr->second);
if (edep > eThreshold) nbOfFired++;
///G4int copyNb = (itr->first);
///G4cout << G4endl << " cryst" << copyNb << ": " << edep/keV << " keV ";
}
if (nbOfFired == 2) fGoodEvents++;
//Dose deposit in patient
//
G4double dose = 0.;
evtMap = static_cast<G4THitsMap<G4double>*>(HCE->GetHC(fCollID_patient));
for (itr = evtMap->GetMap()->begin(); itr != evtMap->GetMap()->end(); itr++) {
///G4int copyNb = (itr->first);
dose = *(itr->second);
}
fSumDose += dose;
fStatDose += dose;
G4Run::RecordEvent(event);
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void Run::Merge(const G4Run* aRun)
{
const Run* localRun = static_cast<const Run*>(aRun);
fGoodEvents += localRun->fGoodEvents;
fSumDose += localRun->fSumDose;
fStatDose += localRun->fStatDose;
G4Run::Merge(aRun);
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
}
Concrete classes of G4VPrimitiveScorer
¶
With Geant4 version 8.0, several concrete primitive scorer classes are
provided, all of which are derived from the G4VPrimitiveScorer
abstract base class and which are to be registered to
G4MultiFunctionalDetector
. Each of them contains one G4THitsMap
object and scores a simple double value for each key.
Track length scorers¶
G4PSTrackLength
The track length is defined as the sum of step lengths of the particles inside the cell. By default, the track weight is not taken into account, but could be used as a multiplier of each step length if the
Weighted()
method of this class object is invoked.G4PSPassageTrackLength
The passage track length is the same as the track length in
G4PSTrackLength
, except that only tracks which pass through the volume are taken into account. It means newly-generated or stopped tracks inside the cell are excluded from the calculation. By default, the track weight is not taken into account, but could be used as a multiplier of each step length if theWeighted()
method of this class object is invoked.
Deposited energy scorers¶
G4PSEnergyDeposit
This scorer stores a sum of particles’ energy deposits at each step in the cell. The particle weight is multiplied at each step.
G4PSDoseDeposit
In some cases, dose is a more convenient way to evaluate the effect of energy deposit in a cell than simple deposited energy. The dose deposit is defined by the sum of energy deposits at each step in a cell divided by the mass of the cell. The mass is calculated from the density and volume of the cell taken from the methods of
G4VSolid
andG4LogicalVolume
. The particle weight is multiplied at each step.G4PSDoseDeposit3D
In the case of replica or nested geometries it is necessary to determine voxel numbers from within the replica hierarchy. For example if the z-axis is parameterised and y is replica of x then the voxel number needs to be calculated accordingly:
G4PSDoseDeposit3D(("DoseDeposit", fNoVoxelsZ, fNoVoxelsY, fNoVoxelsX, 0, 2, 1);
The last three arguments are optional, however required to determine the depth according to each axis (or replica direction). The class creates instances of G4PSDoseDeposit according to:
i * fNj * fNk + j * fNk + k;
where i, j and k correspond to the 3 arguments for number of voxels and replica depth in the declaration.
Current and flux scorers¶
There are two different definitions of a particle’s flow for a given geometry. One is a current and the other is a flux. In our scorers, the current is simply defined as the number of particles (with the particle’s weight) at a certain surface or volume, while the flux takes the particle’s injection angle to the geometry into account. The current and flux are usually defined at a surface, but volume current and volume flux are also provided.
G4PSFlatSurfaceCurrent
Flat surface current is a surface based scorer. The present implementation is limited to scoring only at the -Z surface of a
G4Box
solid. The quantity is defined by the number of tracks that reach the surface. The user must choose a direction of the particle to be scored. The choices are fCurrent_In, fCurrent_Out, or fCurrent_InOut, one of which must be entered as the second argument of the constructor. Here, fCurrent_In scores incoming particles to the cell, while fCurrent_Out scores only outgoing particles from the cell. fCurrent_InOut scores both directions. The current is multiplied by particle weight and is normalized for a unit area.G4PSSphereSurfaceCurrent
Sphere surface current is a surface based scorer, and similar to the G4PSFlatSurfaceCurrent. The only difference is that the surface is defined at the inner surface of a G4Sphere solid.
G4PSPassageCurrent
Passage current is a volume-based scorer. The current is defined by the number of tracks that pass through the volume. A particle weight is applied at the exit point. A passage current is defined for a volume.
G4PSFlatSurfaceFlux
Flat surface flux is a surface based flux scorer. The surface flux is defined by the number of tracks that reach the surface. The expression of surface flux is given by the sum of W/cos(t)/A, where W, t and A represent particle weight, injection angle of particle with respect to the surface normal, and area of the surface. The user must enter one of the particle directions, fFlux_In, fFlux_Out, or fFlux_InOut in the constructor. Here, fFlux_In scores incoming particles to the cell, while fFlux_Out scores outgoing particles from the cell. fFlux_InOut scores both directions.
G4PSCellFlux
Cell flux is a volume based flux scorer. The cell flux is defined by a track length (L) of the particle inside a volume divided by the volume (V) of this cell. The track length is calculated by a sum of the step lengths in the cell. The expression for cell flux is given by the sum of (W*L)/V, where W is a particle weight, and is multiplied by the track length at each step.
G4PSPassageCellFlux
Passage cell flux is a volume based scorer similar to
G4PSCellFlux
. The only difference is that tracks which pass through a cell are taken into account. It means generated or stopped tracks inside the volume are excluded from the calculation.
Other scorers¶
G4PSMinKinEAtGeneration
This scorer records the minimum kinetic energy of secondary particles at their production point in the volume in an event. This primitive scorer does not integrate the quantity, but records the minimum quantity.
G4PSNofSecondary
This class scores the number of secondary particles generated in the volume. The weight of the secondary track is taken into account.
G4PSNofStep
This class scores the number of steps in the cell. A particle weight is not applied.
G4PSCellCharge
This class scored the total charge of particles which has stopped in the volume.
G4VSDFilter
and its derived classes¶
G4VSDFilter
is an abstract class that represents a track filter to
be associated with G4VSensitiveDetector
or G4VPrimitiveScorer
.
It defines a virtual method
G4bool Accept(const G4Step*)
that should return true if this particular step should be scored by
the G4VSensitiveDetector
or G4VPrimitiveScorer
.
While the user can implement his/her own filter class, Geant4 version 8.0 provides the following concrete filter classes:
- G4SDChargedFilter
All charged particles are accepted.
- G4SDNeutralFilter
All neutral particles are accepted.
- G4SDParticleFilter
Particle species which are registered to this filter object by
Add("particle_name")
are accepted. More than one species can be registered.- G4SDKineticEnergyFilter
A track with kinetic energy greater than or equal to EKmin and smaller than EKmin is accepted. EKmin and EKmax should be defined as arguments of the constructor. The default values of EKmin and EKmax are zero and DBL_MAX.
- G4SDParticleWithEnergyFilter
Combination of
G4SDParticleFilter
andG4SDParticleWithEnergyFilter
.
The use of the G4SDParticleFilter
class is demonstrated in
Listing 54, where filters which accept
gamma, electron, positron and electron/positron are defined.
Multiple sensitive detectors associated to a single logical-volume¶
From Geant4 Version 10.3 it is possible to attach multiple sensitive
detectors to a single geometrical element. This is achieved via the use
of a special proxy class, to which multiple child sensitive detectors
are attached: G4MultiSensitiveDetector
. The kernel still sees a single sensitive detector associated to any
given logical-volume, but the proxy will dispatch the calls from kernel
to all the attached child sensitive detectors.
When using the G4VUserDetectorConstruction::SetSensitiveDetector(...)
utility method the handling of multiple sensitive detectors is done
automatically. Multiple calls to the method passing the same logical
volume will trigger the creation and setup of an instance of
G4MultiSensitiveDetector
.
For more complex use cases it may be necessary to manually instantiate
and setup an instance of G4MultiSensitiveDetector
. For this advanced
use case you can refer to the implementation of the
G4VUserDetectorConstruction::SetSensitiveDetector(G4LogicalVolume* logVol, G4VSensitiveDetector* aSD)
utility method.
void MyDetectorConstruction::ConstructSDandField()
{
auto sdman = G4SDManager::GetSDMpointer();
//...
auto mySD = new mySD1("/SD1");
sdman->AddNewDetector(mySD);//Note we explicitly add the SD to the manager
SetSensitiveDetector("LogVolName",mySD);
auto mySD2 = new MySD2("/SD2");
sdman->AddNewDetector(mySD2);
//This will trigger automatic creation and setup of proxy
SetSensitiveDetector("LogVolName",mySD2);
//...
}
Utilities¶
UI-command base scoring¶
Command-based scoring functionality offers the user a possibility to define a scoring-mesh and various scorers for commonly-used physics quantities such as dose, flux, etc. via UI commands.
Due to small performance overhead, it does not come by default. To use this functionality,
G4ScoringManager
has to be activated after the instantiation of G4RunManager in the main()
function, see Listing 57. This will create the UI commands in /score directory.
#include "G4ScoringManager.hh"
int main()
{
// ...
G4RunManager* runManager = new G4RunManager;
G4ScoringManager::GetScoringManager();
// ...
}
An example of a macro creating a scoring mesh of box type with two scorers and a filter is given below:
# Define scoring mesh
/score/create/boxMesh boxMesh_1
/score/mesh/boxSize 100. 100. 100. cm
/score/mesh/nBin 30 30 30
# Define scoring quantity
/score/quantity/energyDeposit boxMash keV
# Define a filter
/score/filter/charged
# Close mesh
/score/close
Detailed usage of command-based scoring is given in the section Command-based scoring.
G4ScoringManager
provides also a default score writer which dumps every
entry of one quantity of a mesh for all quantities of the mesh
one by one in CSV format. To alternate the file format the user
can implement his/her own score writer deriving from G4VUserScoreWriter
base class and set it to G4ScoringManager
. To demonstrate this,
RE03UserScoreWriter
is included in the extended RE03 example in the
runAndEvent category.
Score Ntuple Writer¶
It is also possible to save the scorers hits using Geant4 analysis tools. This functionality is assured by the G4VScoreNtupleWriter
interface (since 10.5) and the G4TScoreNtupleWriter
and G4TScoreNtupleWriterMessenger
classses (since 10.6).
This feature is demonstrated in the basic examples B3 and B4d. The example of activating the score ntuple writer is given below:
#include "G4AnalysisManager.hh"
#include "G4TScoreNtupleWriter.hh"
// Activate score ntuple writer
G4TScoreNtupleWriter<G4AnalysisManager> scoreNtupleWriter;
scoreNtupleWriter.SetVerboseLevel(1);
The Geant4 UI commands defined in G4TScoreNtupleWriterMessenger
can be used to choose
the output file name and the level of verbosity:
/score/ntuple/writerFileName name
/score/ntuple/writerVerbose 1