Belle II Software development
CDCSensitiveDetector.h
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#pragma once
10
11#include <simulation/kernel/SensitiveDetectorBase.h>
12
13#include <cdc/dataobjects/CDCSimHit.h>
14#include <cdc/dataobjects/WireID.h>
15
16#include <mdst/dataobjects/MCParticle.h>
17#include <framework/datastore/StoreArray.h>
18
19#include <vector>
20#include <map>
21
22
23namespace Belle2 {
28 class CDCSimHit;
29 namespace CDC {
30
31 class CDCGeometryPar;
32
34
39
40 public:
41
43 CDCSensitiveDetector(G4String name, G4double thresholdEnergyDeposit, G4double thresholdKineticEnergy);
44
47
49 void Initialize(G4HCofThisEvent*) override;
50
52 bool step(G4Step* aStep, G4TouchableHistory* history) override;
53
55 // void BeginOfEvent(G4HCofThisEvent*);
56
58 void EndOfEvent(G4HCofThisEvent*) override;
59
61 void saveSimHit(const G4int layerId,
62 const G4int wireId,
63 const G4int trackID,
64 const G4int pid,
65 const G4double distance,
66 const G4double tof,
67 const G4double edep,
68 const G4double stepLength,
69 const G4ThreeVector& mom,
70 const G4ThreeVector& posW,
71 const G4ThreeVector& posIn,
72 const G4ThreeVector& posOut,
73 const G4ThreeVector& posTrack,
74 const G4int lr,
75 const G4int NewLrRaw,
76 const G4int NewLr,
77 const G4double speed,
78 const G4double hitWeight);
79
80 //void AddbgOne(bool doit);
81
82 private:
83
85 G4double ClosestApproach(G4ThreeVector bwp, G4ThreeVector fwp, G4ThreeVector posIn, G4ThreeVector posOut,
86 G4ThreeVector& hitPosition, G4ThreeVector& wirePosition);//,G4double& transferT);
87
88
101 void CellBound(const G4int layerId, const G4int ic1, const G4int ic2,
102 const G4double venter[6], const G4double vexit[6],
103 const G4double s1, const G4double s2,
104 G4double xint[6], G4double& sint, G4int& iflag);
105
106
108
122 void GCUBS(const G4double x, const G4double y, const G4double d1, const G4double d2, G4double a[4]);
123
131 void for_Rotat(const G4double bfld[3]);
132
139 void Rotat(G4double& x, G4double& y,
140 G4double& z,
141 const int mode);
142
146 void Rotat(G4double x[3], const int mode);
147
149
174 void HELWIR(const G4double xwb4,
175 const G4double ywb4,
176 const G4double zwb4,
177 const G4double xwf4,
178 const G4double ywf4,
179 const G4double zwf4,
180 const G4double xp,
181 const G4double yp,
182 const G4double zp,
183 const G4double px,
184 const G4double py,
185 const G4double pz,
186 const G4double B_kG[3],
187 const G4double charge,
188 const G4int ntryMax,
189 G4double& distance,
190 G4double q2[3],
191 G4double q1[3],
192 G4double q3[3],
193 G4int& ntry);
194
196
208 void Mvopr(const G4int ndim, const G4double b[3], const G4double m[3][3],
209 const G4double a[3], G4double c[3], const G4int mode);
210
212 std::vector<int> WireId_in_hit_order(int id0, int id1, int nWires);
213
216
219
225 unsigned short areNeighbors(const WireID& wireId, const WireID& otherWireId) const;
226
235 unsigned short areNeighbors(unsigned short iCLayer, unsigned short iSuperLayer, unsigned short iLayer, unsigned short iWire,
236 const WireID& otherWireId) const;
237
238
241
244
249
255
256 G4double m_brot[3][3];
258 // CDCGeometryPar& m_cdcgp; /**< Reference to CDCGeometryPar object. */
265
270
271 G4bool m_wireSag;
279 std::multimap<unsigned short, CDCSimHit*>
282 std::vector<CDCSimHit*> m_hitWithNegWeight;
284 // std::vector<std::multimap<unsigned short, CDCSimHit*>::iterator> m_posWeightMapItBegin; /**< Iterator showing begin of positive hit map*/
285
286 // std::vector<std::multimap<unsigned short, CDCSimHit*>::iterator> m_posWeightMapItEnd; /**< Iterator showing end of positive hit map*/
287
288 // int m_nPosHits, m_nNegHits;
289
290 const signed short CCW = 1;
291 const signed short CW = -1;
292 const signed short CW_OUT_NEIGHBOR = 1;
293 const signed short CW_NEIGHBOR = 3;
294 const signed short CW_IN_NEIGHBOR = 5;
295 const signed short CCW_IN_NEIGHBOR = 7;
296 const signed short CCW_NEIGHBOR = 9;
297 const signed short CCW_OUT_NEIGHBOR = 11;
298 };
299 }
301} // end of namespace Belle2
The Class for CDC Geometry Parameters.
The Class for CDC Sensitive Detector.
const signed short CW_NEIGHBOR
Constant for clockwise.
G4double m_thresholdEnergyDeposit
Threshold energy deposit to be stored.
const signed short CCW_NEIGHBOR
Constant for counterclockwise.
CDCGeometryPar * m_cdcgp
Pointer to CDCGeometryPar object.
const signed short CW_IN_NEIGHBOR
Constant for clockwise inwards.
std::multimap< unsigned short, CDCSimHit * > m_hitWithPosWeight
Map containing hits with positive weight.
const signed short CCW_OUT_NEIGHBOR
Constant for counterclockwise outwards.
G4bool m_wireSag
Switch to activate wire sag effect.
StoreArray< CDCSimHit > m_CDCSimHits
CDC simulation hits.
G4double m_thresholdKineticEnergy
Threshold kinetic energy to be stored.
const signed short CW_OUT_NEIGHBOR
Constant for clockwise outwards.
const signed short CCW
Constant for counterclockwise orientation.
const signed short CCW_IN_NEIGHBOR
Constant for counterclockwise inwards.
G4bool m_modifiedLeftRightFlag
Switch for left/right flag modified for tracking.
G4int m_magneticField
Magnetic field is on or off.
int m_hitNumber
The current number of created hits in an event.
G4double m_brot[3][3]
a rotation matrix.
std::vector< CDCSimHit * > m_hitWithNegWeight
Vector containing hits with negative weight.
G4int m_nonUniformField
Magnetic field is uniform or non-uniform.
StoreArray< MCParticle > m_MCParticles
MC particles.
const signed short CW
Constant for clockwise orientation.
Base class for all Sensitive Detectors to create hits during simulation.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Class to identify a wire inside the CDC.
Definition: WireID.h:34
bool step(G4Step *aStep, G4TouchableHistory *history) override
Process each step and calculate variables defined in CDCB4VHit.
void setModifiedLeftRightFlag()
set left/right flag modified for tracking
unsigned short areNeighbors(const WireID &wireId, const WireID &otherWireId) const
Check if neighboring cell in the same super-layer; essentially a copy from cdcLocalTracking/mclookup.
void EndOfEvent(G4HCofThisEvent *) override
Do what you want to do at the beginning of each event (why this is not called ?)
void for_Rotat(const G4double bfld[3])
Calculates a rotation matrix.
void HELWIR(const G4double xwb4, const G4double ywb4, const G4double zwb4, const G4double xwf4, const G4double ywf4, const G4double zwf4, const G4double xp, const G4double yp, const G4double zp, const G4double px, const G4double py, const G4double pz, const G4double B_kG[3], const G4double charge, const G4int ntryMax, G4double &distance, G4double q2[3], G4double q1[3], G4double q3[3], G4int &ntry)
Calculate closest points between helix and wire.
void CellBound(const G4int layerId, const G4int ic1, const G4int ic2, const G4double venter[6], const G4double vexit[6], const G4double s1, const G4double s2, G4double xint[6], G4double &sint, G4int &iflag)
Calculate intersection of track with cell boundary.
void Rotat(G4double &x, G4double &y, G4double &z, const int mode)
Translation method.
void saveSimHit(const G4int layerId, const G4int wireId, const G4int trackID, const G4int pid, const G4double distance, const G4double tof, const G4double edep, const G4double stepLength, const G4ThreeVector &mom, const G4ThreeVector &posW, const G4ThreeVector &posIn, const G4ThreeVector &posOut, const G4ThreeVector &posTrack, const G4int lr, const G4int NewLrRaw, const G4int NewLr, const G4double speed, const G4double hitWeight)
Save CDCSimHit into datastore.
void Mvopr(const G4int ndim, const G4double b[3], const G4double m[3][3], const G4double a[3], G4double c[3], const G4int mode)
Calculate the result of a matrix times vector.
void Initialize(G4HCofThisEvent *) override
Register CDC hits collection into G4HCofThisEvent.
G4double ClosestApproach(G4ThreeVector bwp, G4ThreeVector fwp, G4ThreeVector posIn, G4ThreeVector posOut, G4ThreeVector &hitPosition, G4ThreeVector &wirePosition)
Assume line track to calculate distance between track and wire (drift length).
void GCUBS(const G4double x, const G4double y, const G4double d1, const G4double d2, G4double a[4])
void reAssignLeftRightInfo()
Re-assign left/right info.
std::vector< int > WireId_in_hit_order(int id0, int id1, int nWires)
Sort wire id.
Abstract base class for different kinds of events.