ATLAS Offline Software
Loading...
Searching...
No Matches
DMCalibrationCalculator.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3*/
4
5// LArG4::DM:CalibrationCalculator
6// Prepared 11-apr-2006 G.Unal
7
8// This class calculates the values needed for calibration hits in the
9// simulation for steps in the Dead Matter between Barrel and End-Cap cryostats
10
11// #define DEBUG_HITS
12
14
16
17#include "G4Step.hh"
18#include "G4StepPoint.hh"
19#include "G4VPhysicalVolume.hh"
20#include "G4ThreeVector.hh"
21#include "globals.hh"
22
23#include <map>
24#include <algorithm>
25#include <cmath>
26#include <climits>
27
28namespace LArG4 {
29
30 namespace DM {
31
33 // Methods
35
36 CalibrationCalculator::CalibrationCalculator(const std::string& name, ISvcLocator *pSvcLocator)
37 : LArCalibCalculatorSvcImp(name,pSvcLocator)
38 {
39 }
40
41
45
46 G4bool CalibrationCalculator::Process(const G4Step* step, LArG4Identifier & identifier,
47 std::vector<G4double> & energies,
49 {
50 // Use the calculators to determine the energies and the
51 // identifier associated with this G4Step. Note that the
52 // default is to process both the energy and the ID.
53
55 {
56 m_energyCalculator.Energies( step, energies );
57 }
58 else
59 for (unsigned int i=0; i != 4; i++) energies.push_back( 0. );
60
61
62 identifier.clear();
63 if ( process == kEnergyAndID || process == kOnlyID )
64 {
65 // Calculate the identifier.
66
67 // static ?
68 static const double oneOverDeta = 10.; // 1/Deta = 1./0.1 = 10.
69 static const double oneOverDphi = 32./M_PI; // 1/Dphi
70
71 // Calculate the mid-point of the step, and the simple geometry variables.
72
73 G4StepPoint* pre_step_point = step->GetPreStepPoint();
74 G4StepPoint* post_step_point = step->GetPostStepPoint();
75
76 G4ThreeVector startPoint = pre_step_point->GetPosition();
77 G4ThreeVector endPoint = post_step_point->GetPosition();
78 G4ThreeVector p = (startPoint + endPoint) * 0.5;
79
80 //G4double rho = p.perp();
81 G4double eta = fabs( p.pseudoRapidity() );
82 G4double phi = p.phi();
83 if ( phi < 0. ) phi += 2.*M_PI; // Normalize for phiBin calculation
84
85#ifdef DEBUG_HITS
86 ATH_MSG_DEBUG("Process: rho,eta,phi " << p.perp() << " " << eta << " " << phi);
87#endif
88
89 // subdet = +/-4 "+" or " -" according to sign of Z in World coorinate
90 G4int subdet = ( p.z() > 0.) ? 4 : -4;
91 G4int phiBin = (int) ( phi * oneOverDphi );
92 if (phiBin>63) phiBin=63;
93 G4int type = 1;
94
95 // Initialize identifier variables with (invalid) default
96 // values (INT_MIN is defined in <climits>).
97 G4int sampling = INT_MIN;
98 G4int region = INT_MIN;
99 G4int etaBin = INT_MIN;
100
101 if (eta < 1.0) {
102 sampling = 2;
103 region = 1;
104 etaBin = (int) ( eta * oneOverDeta );
105 }
106 else if (eta < 1.5) {
107 sampling = 2;
108 region = 2;
109 etaBin = (int) ( (eta-1.) * oneOverDeta );
110 }
111 else if (eta < 1.6) {
112 sampling = 1;
113 region = 4;
114 etaBin = 0;
115 }
116 else if (eta < 1.8) {
117 sampling = 1;
118 region = 5;
119 etaBin = (int) ( (eta-1.5) * oneOverDeta );
120 }
121 else if (eta < 3.2) {
122 sampling = 1;
123 region = 6;
124 etaBin = (int) ( (eta-1.3) * oneOverDeta );
125 }
126 else {
127 ATH_MSG_WARNING("hit above 3.2 in eta !!!");
128 return false;
129 }
130
131
132 // Append the cell ID to the (empty) identifier.
133 identifier << 10 // Calorimeter
134 << subdet // LAr +/-4 where "+" or " -" according to
135 // the sign of Z in World coorinate
136 << type
137 << sampling
138 << region
139 << etaBin
140 << phiBin;
141
142 }
143
144#ifdef DEBUG_HITS
145 G4double energy = accumulate(energies.begin(),energies.end(),0.);
146 ATH_MSG_DEBUG("Process:"
147 << " ID=" << std::string(identifier)
148 << " energy=" << energy
149 << " energies=(" << energies[0]
150 << "," << energies[1]
151 << "," << energies[2]
152 << "," << energies[3] << ")");
153#endif
154
155 // Check for bad result.
156 if ( identifier == LArG4Identifier() )
157 return false;
158
159 return true;
160 }
161
162 } // namespace DM
163
164} // namespace LArG4
#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.
LArCalibCalculatorSvcImp(const std::string &name, ISvcLocator *pSvcLocator)
CaloG4::SimulationEnergies m_energyCalculator
CalibrationCalculator(const std::string &name, ISvcLocator *pSvcLocator)
virtual G4bool Process(const G4Step *step, LArG4Identifier &identifier, std::vector< G4double > &energies, const eCalculatorProcessing process=kEnergyAndID) const override final
const std::string process
eCalculatorProcessing
@ kOnlyEnergy
@ kEnergyAndID