ATLAS Offline Software
Loading...
Searching...
No Matches
LArG4::DM::CalibrationCalculator Class Reference

#include <DMCalibrationCalculator.h>

Inheritance diagram for LArG4::DM::CalibrationCalculator:
Collaboration diagram for LArG4::DM::CalibrationCalculator:

Public Member Functions

 CalibrationCalculator (const std::string &name, ISvcLocator *pSvcLocator)
virtual ~CalibrationCalculator ()
virtual G4bool Process (const G4Step *step, LArG4Identifier &identifier, LArG4Identifier &identifier_sr, std::vector< double > &energies, const LArG4::eCalculatorProcessing process) const override final

Private Attributes

CaloG4::SimulationEnergies m_energyCalculator

Detailed Description

Definition at line 34 of file DMCalibrationCalculator.h.

Constructor & Destructor Documentation

◆ CalibrationCalculator()

LArG4::DM::CalibrationCalculator::CalibrationCalculator ( const std::string & name,
ISvcLocator * pSvcLocator )

Definition at line 36 of file DMCalibrationCalculator.cxx.

37 : LArCalibCalculatorSvcImp(name,pSvcLocator)
38 {
39 }
LArCalibCalculatorSvcImp(const std::string &name, ISvcLocator *pSvcLocator)

◆ ~CalibrationCalculator()

LArG4::DM::CalibrationCalculator::~CalibrationCalculator ( )
virtual

Definition at line 42 of file DMCalibrationCalculator.cxx.

43 {
44 }

Member Function Documentation

◆ Process()

G4bool LArG4::DM::CalibrationCalculator::Process ( const G4Step * step,
LArG4Identifier & identifier,
LArG4Identifier & identifier_sr,
std::vector< double > & energies,
const LArG4::eCalculatorProcessing process ) const
finaloverridevirtual

Definition at line 46 of file DMCalibrationCalculator.cxx.

51 {
52 // Use the calculators to determine the energies and the
53 // identifier associated with this G4Step. Note that the
54 // default is to process both the energy and the ID.
55
57 {
58 m_energyCalculator.Energies( step, energies );
59 }
60 else
61 for (unsigned int i=0; i != 4; i++) energies.push_back( 0. );
62
63
64 identifier.clear();
65 if ( process == kEnergyAndID || process == kOnlyID )
66 {
67 // Calculate the identifier.
68
69 // static ?
70 static const double oneOverDeta = 10.; // 1/Deta = 1./0.1 = 10.
71 static const double oneOverDphi = 32./M_PI; // 1/Dphi
72
73 // Calculate the mid-point of the step, and the simple geometry variables.
74
75 G4StepPoint* pre_step_point = step->GetPreStepPoint();
76 G4StepPoint* post_step_point = step->GetPostStepPoint();
77
78 G4ThreeVector startPoint = pre_step_point->GetPosition();
79 G4ThreeVector endPoint = post_step_point->GetPosition();
80 G4ThreeVector p = (startPoint + endPoint) * 0.5;
81
82 //G4double rho = p.perp();
83 G4double eta = fabs( p.pseudoRapidity() );
84 G4double phi = p.phi();
85 if ( phi < 0. ) phi += 2.*M_PI; // Normalize for phiBin calculation
86
87#ifdef DEBUG_HITS
88 ATH_MSG_DEBUG("Process: rho,eta,phi " << p.perp() << " " << eta << " " << phi);
89#endif
90
91 // subdet = +/-4 "+" or " -" according to sign of Z in World coorinate
92 G4int subdet = ( p.z() > 0.) ? 4 : -4;
93 G4int phiBin = (int) ( phi * oneOverDphi );
94 if (phiBin>63) phiBin=63;
95 G4int type = 1;
96
97 // Initialize identifier variables with (invalid) default
98 // values (INT_MIN is defined in <climits>).
99 G4int sampling = INT_MIN;
100 G4int region = INT_MIN;
101 G4int etaBin = INT_MIN;
102
103 if (eta < 1.0) {
104 sampling = 2;
105 region = 1;
106 etaBin = (int) ( eta * oneOverDeta );
107 }
108 else if (eta < 1.5) {
109 sampling = 2;
110 region = 2;
111 etaBin = (int) ( (eta-1.) * oneOverDeta );
112 }
113 else if (eta < 1.6) {
114 sampling = 1;
115 region = 4;
116 etaBin = 0;
117 }
118 else if (eta < 1.8) {
119 sampling = 1;
120 region = 5;
121 etaBin = (int) ( (eta-1.5) * oneOverDeta );
122 }
123 else if (eta < 3.2) {
124 sampling = 1;
125 region = 6;
126 etaBin = (int) ( (eta-1.3) * oneOverDeta );
127 }
128 else {
129 ATH_MSG_WARNING("hit above 3.2 in eta !!!");
130 return false;
131 }
132
133
134 // Append the cell ID to the (empty) identifier.
135 identifier << 10 // Calorimeter
136 << subdet // LAr +/-4 where "+" or " -" according to
137 // the sign of Z in World coorinate
138 << type
139 << sampling
140 << region
141 << etaBin
142 << phiBin;
143
144 }
145
146#ifdef DEBUG_HITS
147 G4double energy = accumulate(energies.begin(),energies.end(),0.);
148 ATH_MSG_DEBUG("Process:"
149 << " ID=" << std::string(identifier)
150 << " energy=" << energy
151 << " energies=(" << energies[0]
152 << "," << energies[1]
153 << "," << energies[2]
154 << "," << energies[3] << ")");
155#endif
156
157 // Check for bad result.
158 if ( identifier == LArG4Identifier() )
159 return false;
160
161 return true;
162 }
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
bool accumulate(AccumulateMap &map, std::vector< module_t > const &modules, FPGATrackSimMatrixAccumulator const &acc)
Accumulates an accumulator (e.g.
CaloG4::SimulationEnergies m_energyCalculator
const std::string process
@ kOnlyEnergy
@ kEnergyAndID
setSAddress setEtaMS setDirPhiMS setDirZMS setBarrelRadius setEndcapAlpha setEndcapRadius setPhiMap phiBin
setSAddress setEtaMS setDirPhiMS setDirZMS setBarrelRadius setEndcapAlpha setEndcapRadius setInterceptInner setEtaMap etaBin

Member Data Documentation

◆ m_energyCalculator

CaloG4::SimulationEnergies LArG4::DM::CalibrationCalculator::m_energyCalculator
private

Definition at line 58 of file DMCalibrationCalculator.h.


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