Belle II Software  release-08-01-10
SensitiveDetector Class Reference

Sensitive Detector implementation of the CSI detector. More...

#include <SensitiveDetector.h>

Inheritance diagram for SensitiveDetector:
Collaboration diagram for SensitiveDetector:

Public Member Functions

 SensitiveDetector ()
 Constructor.
 
 ~SensitiveDetector ()
 Destructor.
 
int saveSimHit (const G4int cellId, const G4int trackID, const G4int pid, const G4double tof, const G4double edep, G4ThreeVector mom, G4ThreeVector pos, G4ThreeVector PosCell, G4ThreeVector VecCell)
 Save ECLSimHits of the event into datastore. More...
 

Static Public Member Functions

static const std::map< std::string, RelationArray::EConsolidationAction > & getMCParticleRelations ()
 Return a list of all registered Relations with MCParticles.
 
static void setActive (bool activeStatus)
 Enable/Disable all Sensitive Detectors. More...
 
static void registerMCParticleRelation (const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
 Register an relation involving MCParticles. More...
 
static void registerMCParticleRelation (const RelationArray &relation, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
 Overload to make it easer to register MCParticle relations. More...
 

Protected Member Functions

bool step (G4Step *aStep, G4TouchableHistory *) override
 Step processing method. More...
 

Private Member Functions

virtual bool ProcessHits (G4Step *aStep, G4TouchableHistory *aROhist)
 Check if recording hits is enabled and if so call step() and set the correct MCParticle flag. More...
 

Private Attributes

int m_hitNum
 members of SensitiveDetector

 
int m_EvnetNumber
 The current number of created hits in an event. More...
 
int m_oldEvnetNumber
 The current number of created hits in an event. More...
 
int m_trackID
 The current number of created hits in an event. More...
 
G4ThreeVector m_startPos
 track id
 
G4ThreeVector m_endPos
 Position of prestep.
 
G4ThreeVector m_WightedPos
 Position of poststep.
 
G4ThreeVector m_momentum
 Wighted step Position.
 
double m_startTime
 momentum of track
 
double m_endTime
 global time
 
double m_WightedTime
 global time
 
double m_startEnergy
 global time
 
double m_energyDeposit
 particle energy at the entrance in volume
 
double m_trackLength
 energy deposited in volume
 
int CsiHitIndex [18][80]
 length of the track in the volume
 
int iECLCell
 Hit index of StoreArray.
 
int TimeIndex
 Hit Energy of StoreArray.
 
double local_pos
 Hit Time of StoreArray.
 
double T_ave
 position alongthe vector of crystal axis

 
int firstcall
 flight time to diode sensor

 
int m_phiID
 flag of first call

 
int m_thetaID
 The current phi ID in an event. More...
 
int m_cellID
 The current theta ID in an event. More...
 
Const::EDetector m_subdetector
 Subdetector the class belongs to.
 

Static Private Attributes

static std::map< std::string, RelationArray::EConsolidationActions_mcRelations
 Static set holding all relations which have to be updated at the end of the Event.
 
static bool s_active
 Static bool which indicates wether recording of hits is enabled.
 

Detailed Description

Sensitive Detector implementation of the CSI detector.

Definition at line 22 of file SensitiveDetector.h.

Member Function Documentation

◆ ProcessHits()

bool ProcessHits ( G4Step *  aStep,
G4TouchableHistory *  aROhist 
)
inlineprivatevirtualinherited

Check if recording hits is enabled and if so call step() and set the correct MCParticle flag.

Called by Geant4 for each step inside the sensitive volumes attached

Definition at line 94 of file SensitiveDetectorBase.h.

◆ registerMCParticleRelation() [1/2]

static void registerMCParticleRelation ( const RelationArray relation,
RelationArray::EConsolidationAction  ignoreAction = RelationArray::c_negativeWeight 
)
inlinestaticinherited

Overload to make it easer to register MCParticle relations.

Parameters
relationRelationArray to register
ignoreAction

Definition at line 66 of file SensitiveDetectorBase.h.

◆ registerMCParticleRelation() [2/2]

void registerMCParticleRelation ( const std::string &  name,
RelationArray::EConsolidationAction  ignoreAction = RelationArray::c_negativeWeight 
)
staticinherited

Register an relation involving MCParticles.

All Relations which point from an MCParticle to something have to be registered with addMCParticleRelation() because the index of the MCParticles might change at the end of the event. During simulation, the TrackID should be used as index of the MCParticle

Parameters
nameName of the relation to register
ignoreAction

Definition at line 22 of file SensitiveDetectorBase.cc.

◆ saveSimHit()

int saveSimHit ( const G4int  cellId,
const G4int  trackID,
const G4int  pid,
const G4double  tof,
const G4double  edep,
G4ThreeVector  mom,
G4ThreeVector  pos,
G4ThreeVector  PosCell,
G4ThreeVector  VecCell 
)

Save ECLSimHits of the event into datastore.

Parameters
cellIdthe channel identifier
trackIDthe track identifier
pidthe pdg code of the particle
tofthe time of flight of the particle
edepthe energy deposited in the crystal by the particle
momthe momentum of the particle
posTO BE COMPLETED
PosCellTO BE COMPLETED
VecCellTO BE COMPLETED
Returns
int the number of hits saved in this event

Definition at line 153 of file SensitiveDetector.cc.

163  {
164 
165  //Get the datastore arrays
166  StoreArray<MCParticle> mcParticles;
167  StoreArray<CsiSimHit> simHits;
168 
169  RelationArray simHitRel(mcParticles, simHits);
170 
171  // convert Hep3Vectors to TVector3
172  TVector3 momentum(mom.getX(), mom.getY(), mom.getZ());
173  TVector3 position(pos.getX(), pos.getY(), pos.getZ());
174 
175  //Append
176  simHits.appendNew(cellId, trackID, pid, tof / CLHEP::ns, edep / CLHEP::GeV,
177  momentum * (1.0 / CLHEP::GeV), position * (1.0 / CLHEP::cm));
178 
179  int simhitNumber = simHits.getEntries() - 1;
180  simHitRel.add(trackID, simhitNumber);
181  B2DEBUG(150, "Saving CsiSimHit Number: " << simhitNumber);
182 
183  // Store CsiHit: total energy deposited per hit, with average time.
184 
185  // Old names:
186  /*
187  StoreArray<CsiHit> eclHitArray;
188  RelationArray eclHitRel(mcParticles, eclHitArray);
189  */
190 
191  // New names:
192 
193  StoreArray<CsiHit> csiHits;
194 
195  RelationArray csiHitRel(mcParticles, csiHits);
196 
197 
198  StoreObjPtr<EventMetaData> eventMetaDataPtr;
199  int m_currentEvnetNumber = eventMetaDataPtr->getEvent();
200 
201 
202  if (firstcall == 0 || m_currentEvnetNumber != m_oldEvnetNumber) {
203  m_oldEvnetNumber = m_currentEvnetNumber;
204  for (int iCSICell = 0; iCSICell < 18; iCSICell++) {
205  for (int TimeInd = 0; TimeInd < 80; TimeInd++) {
206  CsiHitIndex[iCSICell][TimeInd] = -1;
207  }
208  }
209  firstcall++;
210  }
211 
212 
213  if (m_currentEvnetNumber == m_oldEvnetNumber) {
214  if ((tof / CLHEP::ns) < 8000) {
215  TimeIndex = (int)((tof / CLHEP::ns) / 100);
216 
217  // local_pos is the distance between the hit and the PIN-diode end of the crystal.
218  local_pos = (15. - (pos - PosCell) * VecCell);
219 
220  // From the ECL package T_ave = 6.05 + 0.0749 * local_pos - 0.00112 * local_pos * local_pos + (tof / CLHEP::ns) ;
221 
222  //For CsI: we remove the offset (I hypothetize it if from photo-sensor delay)
223  // we remove the acceleration of the light
224  // we correct the velocity term using n_CsI=1.80
225  T_ave = 0.0600 * local_pos + (tof / CLHEP::ns) ;
226 
227 
228  double E_cell = (edep / CLHEP::GeV);
229  if (CsiHitIndex[cellId][TimeIndex] == -1) {
230  //m_hitNum = csiHits->GetLast() + 1;
231  //new(csiHits->AddrAt(m_hitNum)) CSIHit();
232  csiHits.appendNew();
233  m_hitNum = csiHits.getEntries() - 1;
234  CsiHitIndex[cellId][TimeIndex] = m_hitNum;
235  csiHits[m_hitNum]->setCellId(cellId);
236  csiHits[m_hitNum]->setEnergyDep(E_cell);
237  csiHits[m_hitNum]->setTimeAve(T_ave);
238  csiHits[m_hitNum]->setTimeVar(0);
239  csiHitRel.add(trackID, m_hitNum);
240 
241  } else {
242  m_hitNum = CsiHitIndex[cellId][TimeIndex];
243  double old_edep = csiHits[m_hitNum]->getEnergyDep();
244  double old_TimeAve = csiHits[m_hitNum]->getTimeAve();
245  double old_TimeVar = csiHits[m_hitNum]->getTimeVar();
246 
247  //double new_edep = E_cell + old_edep;
248  //double new_TimeAve = (old_edep * old_TimeAve + E_cell * T_ave) / new_edep;
249 
250  double temp = E_cell + old_edep;
251  double delta = T_ave - old_TimeAve;
252  double R = delta * E_cell / temp;
253  double new_TimeAve = old_TimeAve + R;
254  double new_TimeVar = old_TimeVar + old_edep * delta * R;
255  double new_edep = temp;
256 
257  //Algo from Wiki:
258  /*
259  for x, weight in dataWeightPairs: # Alternatively "for x, weight in zip(data, weights):"
260  temp = weight + sumweight
261  delta = x - mean
262  R = delta * weight / temp
263  mean = mean + R
264  M2 = M2 + sumweight * delta * R # Alternatively, "M2 = M2 + weight * delta * (x−mean)"
265  sumweight = temp
266 
267  variance_n = M2/sumweight
268  */
269 
270  csiHits[m_hitNum]->setEnergyDep(new_edep);
271  csiHits[m_hitNum]->setTimeAve(new_TimeAve);
272  csiHits[m_hitNum]->setTimeVar(new_TimeVar);
273  csiHitRel.add(trackID, m_hitNum);
274 
275  }
276  }
277  } else { B2ERROR("m_currentEvnetNumber ERROR: m_oldEvnetNumber==m_oldEvnetNumber"); }
278 
279  return (simhitNumber);
280  }//saveSimHit
double R
typedef autogenerated by FFTW
int TimeIndex
Hit Energy of StoreArray.
int m_hitNum
members of SensitiveDetector
int firstcall
flight time to diode sensor
double T_ave
position alongthe vector of crystal axis
double local_pos
Hit Time of StoreArray.
int CsiHitIndex[18][80]
length of the track in the volume
int m_oldEvnetNumber
The current number of created hits in an event.

◆ setActive()

static void setActive ( bool  activeStatus)
inlinestaticinherited

Enable/Disable all Sensitive Detectors.

By default, all sensitive detectors won't create hits to make it possible to use the Geant4 Navigator for non-simulation purposes. Only during simulation the sensitive detectors will be enabled to record hits

Parameters
activeStatusbool to indicate wether hits should be recorded

Definition at line 50 of file SensitiveDetectorBase.h.

◆ step()

bool step ( G4Step *  aStep,
G4TouchableHistory *   
)
overrideprotectedvirtual

Step processing method.

Parameters
aStepthe G4Step with the current step information
Returns
true if a Hit has been created, false if the hit was ignored

Implements SensitiveDetectorBase.

Definition at line 85 of file SensitiveDetector.cc.

Member Data Documentation

◆ m_cellID

int m_cellID
private

The current theta ID in an event.

Used to fill the DataStore ECL array

Definition at line 93 of file SensitiveDetector.h.

◆ m_EvnetNumber

int m_EvnetNumber
private

The current number of created hits in an event.

Used to fill the DataStore ECLHit.

Definition at line 69 of file SensitiveDetector.h.

◆ m_oldEvnetNumber

int m_oldEvnetNumber
private

The current number of created hits in an event.

Used to fill the DataStore ECL EB array.

Definition at line 70 of file SensitiveDetector.h.

◆ m_thetaID

int m_thetaID
private

The current phi ID in an event.

Used to fill the DataStore ECL array

Definition at line 92 of file SensitiveDetector.h.

◆ m_trackID

int m_trackID
private

The current number of created hits in an event.

Used to fill the DataStore

Definition at line 71 of file SensitiveDetector.h.


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