1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
|
#ifndef EventFilter_SiPixelRawToDigi_interface_PixelUnpackingRegions_h
#define EventFilter_SiPixelRawToDigi_interface_PixelUnpackingRegions_h
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/ESWatcher.h"
#include "FWCore/Utilities/interface/ESGetToken.h"
#include "DataFormats/Math/interface/Point3D.h"
#include "DataFormats/Math/interface/Vector3D.h"
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#include "DataFormats/Candidate/interface/LeafCandidate.h"
#include "CondFormats/DataRecord/interface/SiPixelFedCablingMapRcd.h"
#include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingMap.h"
#include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingTree.h"
#include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
#include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
#include "FWCore/Framework/interface/ConsumesCollector.h"
#include <cmath>
#include <vector>
#include <set>
/** \class PixelUnpackingRegions
*
* Input: One or several collections of Candidate-based seeds with their objects
* defining the directions of unpacking regions; separate deltaPhi and maxZ
* tolerances could be given to each input collection.
* Output: FED ids and module detIds that need to be unpacked
*
*/
class PixelUnpackingRegions {
public:
/// container to define regions for objects of interest in each event by:
/// object direction
/// dphi max distance from region direction to center of a pixel module
/// maxZ max projected z of a pixel module (when projecting along region direction onto beamline)
struct Region {
Region(const math::XYZVector& dir, float dphi = 0.5f, float maxz = 24.f) : v(dir), dPhi(dphi), maxZ(maxz) {
cosphi = v.x() / v.rho();
sinphi = v.y() / v.rho();
atantheta = v.z() / v.rho();
}
math::XYZVector v;
float dPhi, maxZ;
float cosphi, sinphi, atantheta;
};
PixelUnpackingRegions(const edm::ParameterSet&, edm::ConsumesCollector&& iC);
~PixelUnpackingRegions() {}
/// has to be run during each event
void run(const edm::Event& e, const edm::EventSetup& es);
/// check whether a FED has to be unpacked
bool mayUnpackFED(unsigned int fed_n) const;
/// check whether a module has to be unpacked
bool mayUnpackModule(unsigned int id) const;
/// full set of module ids to unpack
const std::set<unsigned int>* modulesToUnpack() const { return &modules_; }
/// various informational accessors:
unsigned int nFEDs() const { return feds_.size(); }
unsigned int nModules() const { return modules_.size(); }
unsigned int nBarrelModules() const;
unsigned int nForwardModules() const;
unsigned int nRegions() const { return nreg_; }
struct Module {
float phi;
float x, y, z;
unsigned int id;
unsigned int fed;
Module() {}
Module(float ph) : phi(ph), x(0.f), y(0.f), z(0.f), id(0), fed(0) {}
bool operator<(const Module& m) const {
if (phi < m.phi)
return true;
if (phi == m.phi && id < m.id)
return true;
return false;
}
};
private:
// input parameters
std::vector<edm::InputTag> inputs_;
std::vector<double> dPhi_;
std::vector<double> maxZ_;
edm::InputTag beamSpotTag_;
edm::EDGetTokenT<reco::BeamSpot> tBeamSpot;
std::vector<edm::EDGetTokenT<reco::CandidateView>> tCandidateView;
std::set<unsigned int> feds_;
std::set<unsigned int> modules_;
unsigned int nreg_;
/// run by the run method: (re)initialize the cabling data when it's necessary
void initialize(const edm::EventSetup& es);
// add a new direction of a region of interest
void addRegion(Region& r);
// gather info into feds_ and modules_ from a range of a Module vector
void gatherFromRange(Region& r, std::vector<Module>::const_iterator, std::vector<Module>::const_iterator);
// addRegion for a local (BPIX or +-FPIX) container
void addRegionLocal(Region& r, std::vector<Module>& container, const Module& lo, const Module& hi);
// local containers of barrel and endcaps Modules sorted by phi
std::vector<Module> phiBPIX_;
std::vector<Module> phiFPIXp_;
std::vector<Module> phiFPIXm_;
std::unique_ptr<SiPixelFedCablingTree> cabling_;
math::XYZPoint beamSpot_;
edm::ESWatcher<SiPixelFedCablingMapRcd> watcherSiPixelFedCablingMap_;
edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> trackerGeomToken_;
edm::ESGetToken<SiPixelFedCablingMap, SiPixelFedCablingMapRcd> cablingMapToken_;
};
#endif // EventFilter_SiPixelRawToDigi_interface_PixelUnpackingRegions_h
|