Geant4 User's Guide
For Application Developers Detector Definition and Response |
4.4 Hits
4.4.1 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
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.
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 Section 8.9.
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
Source listing 4.4.1 shows an example of a concrete hit class.
#ifndef ExN04TrackerHit_h #define ExN04TrackerHit_h 1 #include "G4VHit.hh" #include "G4THitsCollection.hh" #include "G4Allocator.hh" #include "G4ThreeVector.hh" class ExN04TrackerHit : public G4VHit { public: ExN04TrackerHit(); ~ExN04TrackerHit(); ExN04TrackerHit(const ExN04TrackerHit &right); const ExN04TrackerHit& operator=(const ExN04TrackerHit &right); int operator==(const ExN04TrackerHit &right) const; inline void * operator new(size_t); inline void operator delete(void *aHit); void Draw() const; void Print() const; private: G4double edep; G4ThreeVector pos; public: inline void SetEdep(G4double de) { edep = de; } inline G4double GetEdep() const { return edep; } inline void SetPos(G4ThreeVector xyz) { pos = xyz; } inline G4ThreeVector GetPos() const { return pos; } }; typedef G4THitsCollection<ExN04TrackerHit> ExN04TrackerHitsCollection; extern G4Allocator<ExN04TrackerHit> ExN04TrackerHitAllocator; inline void* ExN04TrackerHit::operator new(size_t) { void *aHit; aHit = (void *) ExN04TrackerHitAllocator.MallocSingle(); return aHit; } inline void ExN04TrackerHit::operator delete(void *aHit) { ExN04TrackerHitAllocator.FreeSingle((ExN04TrackerHit*) aHit); } #endif |
Source listing 4.4.1 An example of a concrete hit class. |
G4Allocator is a class for fast allocation of objects to the heap through the paging mechanism. For details of G4Allocator, refer to 3.2.4. 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 section 4.4.5 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.
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. In the case of a "Readout" geometry (see Section 4.4.3), objects of the G4TouchableHistory class may be used as an optional input.
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 pointer to your sensitive detector must be set to one or more G4LogicalVolume objects to set the sensitivity of these volumes. The pointer should also be registered to G4SDManager, as described in Section 4.4.4.
G4VSensitiveDetector has three major virtual methods.
As an example, the accordion calorimeter of ATLAS has a complicated tracking geometry, however the readout can be done by simple cylindrical sectors divided by theta, phi, and depth. Tracks will be traced in the tracking geometry, the ``real'' one, and the sensitive detector will have its own readout geometry Geant4 will message to find to which ``readout'' cell the current hit belongs.
The picture shows how this association is done in Geant4. The first step is to associate a sensitive detector to a volume of the tracking geometry, in the usual way (see Section 4.4.2). The next step is to associate your G4VReadoutGeometry object to the sensitive detector.
At tracking time, the base class G4VReadoutGeometry will provide to your sensitive detector code the G4TouchableHistory in the Readout geometry at the beginning of the step position (position of PreStepPoint of G4Step) and at this position only.
This G4TouchableHistory is given to your sensitive detector code through the G4VSensitiveDetector virtual method:
G4bool processHits(G4Step* aStep, G4TouchableHistory* ROhist);by the ROhist argument.
Thus, you will be able to use information from both the G4Step and the G4TouchableHistory coming from your Readout geometry. Note that since the association is done through a sensitive detector object, it is perfectly possible to have several Readout geometries in parallel.
Definition of a virtual geometry setup
The base class for the implementation of a Readout geometry is G4VReadoutGeometry. This class has a single pure virtual protected method:
virtual G4VPhysicalVolume* build() = 0;which you must override in your concrete class. The G4VPhysicalVolume pointer you will have to return is of the physical world of the Readout geometry.
The step by step procedure for constructing a Readout geometry is:
The world is specified in the same way as for the detector construction: a physical volume with no mother. The axis system of this world is the same as the one of the world for tracking.
In this geometry you need to declare the sensitive parts in the same way as in the tracking geometry: by setting a non-NULL G4VSensitiveDetector pointer in, say, the relevant G4LogicalVolume objects. This sensitive class needs to be there, but will not be used.
You will also need to assign well defined materials to the volumes you place in this geometry, but these materials are irrelevant since they will not be seen by the tracking. It is foreseen to allow the setting of a NULL pointer in this case of the parallel geometry.
MyROGeom* ROgeom = new MyROGeom("ROName");
ROgeom->buildROGeometry();That will invoke your build() method.
MySensitive->SetROgeometry(ROgeom);
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 myCalwill inactivate all detectors belonging to the myCal category.
Access to the hit collections
Hit collections are accessed for various cases.
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));
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 the next section, 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 next section 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.
Source listing 4.4.2 shows an example of primitive sensitivity class definitions.
void ExN07DetectorConstruction::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 = calName[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); } } } } |
Source listing 4.4.2 An example of defining primitive sensitivity classes taken from ExN07DetectorConstruction. |
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 astarisk ("*"). 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. Source listing 4.4.3 shows the use of G4THitsMap.
#include "ExN07Run.hh" #include "G4Event.hh" #include "G4HCofThisEvent.hh" #include "G4SDManager.hh" ExN07Run::ExN07Run() { G4String detName[6] = {"Calor-A_abs","Calor-A_gap","Calor-B_abs","Calor-B_gap","Calor-C_abs","Calor-C_gap"}; G4String primNameSum[6] = {"eDep","nGamma","nElectron","nPositron","trackLength","nStep"}; G4String primNameMin[3] = {"minEkinGamma","minEkinElectron","minEkinPositron"}; G4SDManager* SDMan = G4SDManager::GetSDMpointer(); G4String fullName; for(size_t i=0;i<6;i++) { for(size_t j=0;j<6;j++) { fullName = detName[i]+"/"+primNameSum[j]; colIDSum[i][j] = SDMan->GetCollectionID(fullName); } for(size_t k=0;k<3;k++) { fullName = detName[i]+"/"+primNameMin[k]; colIDMin[i][k] = SDMan->GetCollectionID(fullName); } } } void ExN07Run::RecordEvent(const G4Event* evt) { G4HCofThisEvent* HCE = evt->GetHCofThisEvent(); if(!HCE) return; numberOfEvent++; for(size_t i=0;i<6;i++) { for(size_t j=0;j<6;j++) { G4THitsMap<G4double>* evtMap = (G4THitsMap<G4double>*)(HCE->GetHC(colIDSum[i][j])); mapSum[i][j] += *evtMap; } for(size_t k=0;k<3;k++) { G4THitsMap<G4double>* evtMap = (G4THitsMap<G4double>*)(HCE->GetHC(colIDMin[i][k])); std::map<G4int,G4double*>::iterator itr = evtMap->GetMap()->begin(); for(; itr != evtMap->GetMap()->end(); itr++) { G4int key = (itr->first); G4double val = *(itr->second); G4double* mapP = mapMin[i][k][key]; if( mapP && (val>*mapP) ) continue; mapMin[i][k].set(key,val); } } } } |
Source listing 4.4.3 An example of accessing to G4THitsMap objects. |
Track length scorers
Deposited energy scorers
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.
Other scorers
While the user can implement his/her own filter class, Geant4 version 8.0 provides the following concrete filter classes:
The use of the G4SDParticleFilter class is demonstrated in the Source
Listing 4.4.2, where filters which accept gamma, electron, positron and
electron/positron are defined.