Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-25 02:13:23

0001 /** \class PixelBaryCentreAnalyzer
0002  *  The analyer works as the following :
0003  *  - Read global tracker position from global tag
0004  *  - Read tracker alignment constants from different ESsource with different labels
0005  *  - Calculate barycentres for different pixel substructures using global tracker position and alignment constants and store them in trees, one for each ESsource label.
0006  *
0007  *  Python script plotBaryCentre_VS_BeamSpot.py under script dir is used to plot barycentres from alignment constants used in Prompt-Reco, End-of-Year Rereco and so-called Run-2 (Ultra)Legacy Rereco. Options of the plotting script can be found from the helper in the script.
0008  *
0009  *  $Date: 2021/01/05 $
0010  *  $Revision: 1.0 $
0011  *  \author Tongguang Cheng - Beihang University <tongguang.cheng@cern.ch>
0012  *
0013 */
0014 
0015 // Framework
0016 #include "FWCore/Framework/interface/Frameworkfwd.h"
0017 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0018 #include "FWCore/Framework/interface/Event.h"
0019 #include "FWCore/Framework/interface/MakerMacros.h"
0020 #include "FWCore/Framework/interface/ESWatcher.h"
0021 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0022 #include "FWCore/Utilities/interface/ESGetToken.h"
0023 
0024 // Phase-1 Pixel
0025 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0026 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0027 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0028 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0029 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0030 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0031 
0032 // pixel quality
0033 #include "CondFormats/SiPixelObjects/interface/SiPixelQuality.h"
0034 #include "CondFormats/DataRecord/interface/SiPixelQualityFromDbRcd.h"
0035 // global postion
0036 #include "CondFormats/Alignment/interface/DetectorGlobalPosition.h"
0037 #include "CondFormats/AlignmentRecord/interface/GlobalPositionRcd.h"
0038 // tracker alignment
0039 #include "CondFormats/AlignmentRecord/interface/TrackerAlignmentRcd.h"
0040 #include "CondFormats/Alignment/interface/Alignments.h"
0041 // beamspot
0042 #include "CondFormats/DataRecord/interface/BeamSpotObjectsRcd.h"
0043 #include "CondFormats/BeamSpotObjects/interface/BeamSpotObjects.h"
0044 
0045 // Point and Vector
0046 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0047 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
0048 
0049 // TFileService
0050 #include "FWCore/ServiceRegistry/interface/Service.h"
0051 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0052 
0053 // ROOT
0054 #include "TTree.h"
0055 #include "TString.h"
0056 
0057 //
0058 // class declaration
0059 //
0060 
0061 class PixelBaryCentreAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0062 public:
0063   explicit PixelBaryCentreAnalyzer(const edm::ParameterSet&);
0064   ~PixelBaryCentreAnalyzer() override = default;
0065 
0066   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0067 
0068   struct SimplePoint {
0069     float x, y, z;
0070     SimplePoint(const GlobalPoint& p) : x(p.x()), y(p.y()), z(p.z()){};
0071     SimplePoint() : x(0), y(0), z(0){};
0072   };
0073   static const unsigned int nPixelLayers = 4;
0074   static const unsigned int nPixelDiscs = 3;
0075 
0076 private:
0077   void beginJob() override;
0078   void analyze(const edm::Event&, const edm::EventSetup&) override;
0079   void endJob() override;
0080 
0081   void initBC();
0082   void initBS();
0083 
0084   const bool usePixelQuality_;
0085   int phase_;
0086 
0087   // ----------member data ---------------------------
0088   edm::ESWatcher<BeamSpotObjectsRcd> watcherBS_;
0089   edm::ESWatcher<TrackerAlignmentRcd> watcherTkAlign_;
0090 
0091   // labels of TkAlign tags
0092   std::vector<std::string> bcLabels_;
0093   // labels of beamspot tags
0094   std::vector<std::string> bsLabels_;
0095 
0096   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> trackerGeometryToken_;
0097   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> trackerTopologyToken_;
0098   const edm::ESGetToken<SiPixelQuality, SiPixelQualityFromDbRcd> siPixelQualityToken_;
0099 
0100   const edm::ESGetToken<Alignments, GlobalPositionRcd> gprToken_;
0101   std::map<std::string, edm::ESGetToken<Alignments, TrackerAlignmentRcd>> tkAlignTokens_;
0102   std::map<std::string, edm::ESGetToken<BeamSpotObjects, BeamSpotObjectsRcd>> bsTokens_;
0103 
0104   // tree content
0105   int run_;
0106   int ls_;
0107 
0108   GlobalPoint BS_;
0109 
0110   GlobalPoint PIX_, BPIX_, FPIX_;
0111   GlobalPoint BPIX_Flipped_, BPIX_NonFlipped_, BPIX_DiffFlippedNonFlipped_;
0112 
0113   GlobalPoint BPIXLayer_[nPixelLayers];
0114   GlobalPoint BPIXLayer_Flipped_[nPixelLayers];
0115   GlobalPoint BPIXLayer_NonFlipped_[nPixelLayers];
0116   GlobalPoint BPIXLayer_DiffFlippedNonFlipped_[nPixelLayers];
0117 
0118   GlobalPoint FPIX_plus_, FPIX_minus_;
0119   GlobalPoint FPIXDisks_plus_[nPixelDiscs];
0120   GlobalPoint FPIXDisks_minus_[nPixelDiscs];
0121 
0122   SimplePoint vBS_;
0123 
0124   SimplePoint vPIX_, vBPIX_, vFPIX_;
0125   SimplePoint vBPIX_Flipped_, vBPIX_NonFlipped_, vBPIX_DiffFlippedNonFlipped_;
0126 
0127   SimplePoint vBPIXLayer_[nPixelLayers];
0128   SimplePoint vBPIXLayer_Flipped_[nPixelLayers];
0129   SimplePoint vBPIXLayer_NonFlipped_[nPixelLayers];
0130   SimplePoint vBPIXLayer_DiffFlippedNonFlipped_[nPixelLayers];
0131 
0132   SimplePoint vFPIX_plus_, vFPIX_minus_;
0133   SimplePoint vFPIXDisks_plus_[nPixelDiscs];
0134   SimplePoint vFPIXDisks_minus_[nPixelDiscs];
0135 
0136   edm::Service<TFileService> tFileService;
0137   std::map<std::string, TTree*> bcTrees_;
0138   std::map<std::string, TTree*> bsTrees_;
0139 };
0140 
0141 //
0142 // constructors and destructor
0143 //
0144 PixelBaryCentreAnalyzer::PixelBaryCentreAnalyzer(const edm::ParameterSet& iConfig)
0145     : usePixelQuality_(iConfig.getUntrackedParameter<bool>("usePixelQuality")),
0146       bcLabels_(iConfig.getUntrackedParameter<std::vector<std::string>>("tkAlignLabels")),
0147       bsLabels_(iConfig.getUntrackedParameter<std::vector<std::string>>("beamSpotLabels")),
0148       trackerGeometryToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()),
0149       trackerTopologyToken_(esConsumes<TrackerTopology, TrackerTopologyRcd>()),
0150       siPixelQualityToken_(esConsumes<SiPixelQuality, SiPixelQualityFromDbRcd>()),
0151       gprToken_(esConsumes<Alignments, GlobalPositionRcd>()) {
0152   for (const auto& label : bcLabels_) {
0153     bcTrees_[label] = nullptr;
0154     tkAlignTokens_[label] = esConsumes<Alignments, TrackerAlignmentRcd>(edm::ESInputTag{"", label});
0155   }
0156 
0157   for (const auto& label : bsLabels_) {
0158     bsTrees_[label] = nullptr;
0159     bsTokens_[label] = esConsumes<BeamSpotObjects, BeamSpotObjectsRcd>(edm::ESInputTag{"", label});
0160   }
0161 
0162   usesResource("TFileService");
0163 }
0164 
0165 //
0166 // member functions
0167 //
0168 
0169 void PixelBaryCentreAnalyzer::initBS() {
0170   double dummy_float = 999999.0;
0171 
0172   BS_ = GlobalPoint(dummy_float, dummy_float, dummy_float);
0173   vBS_ = SimplePoint(BS_);
0174 }
0175 
0176 void PixelBaryCentreAnalyzer::initBC() {
0177   // init to large number (unreasonable number) not zero
0178   double dummy_float = 999999.0;
0179 
0180   PIX_ = GlobalPoint(dummy_float, dummy_float, dummy_float);
0181   BPIX_ = GlobalPoint(dummy_float, dummy_float, dummy_float);
0182   FPIX_ = GlobalPoint(dummy_float, dummy_float, dummy_float);
0183 
0184   BPIX_Flipped_ = GlobalPoint(dummy_float, dummy_float, dummy_float);
0185   BPIX_NonFlipped_ = GlobalPoint(dummy_float, dummy_float, dummy_float);
0186   BPIX_DiffFlippedNonFlipped_ = GlobalPoint(dummy_float, dummy_float, dummy_float);
0187 
0188   FPIX_plus_ = GlobalPoint(dummy_float, dummy_float, dummy_float);
0189   FPIX_minus_ = GlobalPoint(dummy_float, dummy_float, dummy_float);
0190 
0191   for (unsigned int i = 0; i < nPixelLayers; i++) {
0192     BPIXLayer_[i] = GlobalPoint(dummy_float, dummy_float, dummy_float);
0193     BPIXLayer_Flipped_[i] = GlobalPoint(dummy_float, dummy_float, dummy_float);
0194     BPIXLayer_NonFlipped_[i] = GlobalPoint(dummy_float, dummy_float, dummy_float);
0195     BPIXLayer_DiffFlippedNonFlipped_[i] = GlobalPoint(dummy_float, dummy_float, dummy_float);
0196   }
0197 
0198   for (unsigned int i = 0; i < nPixelDiscs; i++) {
0199     FPIXDisks_plus_[i] = GlobalPoint(dummy_float, dummy_float, dummy_float);
0200     FPIXDisks_minus_[i] = GlobalPoint(dummy_float, dummy_float, dummy_float);
0201   }
0202 
0203   vPIX_ = SimplePoint(PIX_);
0204   vBPIX_ = SimplePoint(BPIX_);
0205   vFPIX_ = SimplePoint(FPIX_);
0206 
0207   vBPIX_Flipped_ = SimplePoint(BPIX_Flipped_);
0208   vBPIX_NonFlipped_ = SimplePoint(BPIX_NonFlipped_);
0209   vBPIX_DiffFlippedNonFlipped_ = SimplePoint(BPIX_DiffFlippedNonFlipped_);
0210 
0211   vFPIX_plus_ = SimplePoint(FPIX_plus_);
0212   vFPIX_minus_ = SimplePoint(FPIX_minus_);
0213 
0214   for (unsigned int i = 0; i < nPixelLayers; i++) {
0215     vBPIXLayer_[i] = SimplePoint(BPIXLayer_[i]);
0216     vBPIXLayer_Flipped_[i] = SimplePoint(BPIXLayer_Flipped_[i]);
0217     vBPIXLayer_NonFlipped_[i] = SimplePoint(BPIXLayer_NonFlipped_[i]);
0218     vBPIXLayer_DiffFlippedNonFlipped_[i] = SimplePoint(BPIXLayer_DiffFlippedNonFlipped_[i]);
0219   }
0220 
0221   for (unsigned int i = 0; i < nPixelDiscs; i++) {
0222     vFPIXDisks_plus_[i] = SimplePoint(FPIXDisks_plus_[i]);
0223     vFPIXDisks_minus_[i] = SimplePoint(FPIXDisks_minus_[i]);
0224   }
0225 }
0226 
0227 // ------------ method called for each event  ------------
0228 void PixelBaryCentreAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0229   bool prepareTkAlign = false;
0230   bool prepareBS = false;
0231 
0232   // ES watcher can noly run once in the same event,
0233   // otherwise it will turn false whatsoever because the condition doesn't change in the second time call.
0234   if (watcherTkAlign_.check(iSetup))
0235     prepareTkAlign = true;
0236   if (watcherBS_.check(iSetup))
0237     prepareBS = true;
0238 
0239   if (!prepareTkAlign && !prepareBS)
0240     return;
0241 
0242   run_ = iEvent.id().run();
0243   ls_ = iEvent.id().luminosityBlock();
0244 
0245   if (prepareTkAlign) {  // check for new IOV for TKAlign
0246 
0247     phase_ = -1;
0248 
0249     const TrackerGeometry* tkGeo = &iSetup.getData(trackerGeometryToken_);
0250     const TrackerTopology* tkTopo = &iSetup.getData(trackerTopologyToken_);
0251 
0252     if (tkGeo->isThere(GeomDetEnumerators::PixelBarrel) && tkGeo->isThere(GeomDetEnumerators::PixelEndcap)) {
0253       phase_ = 0;
0254     } else if (tkGeo->isThere(GeomDetEnumerators::P1PXB) && tkGeo->isThere(GeomDetEnumerators::P1PXEC)) {
0255       phase_ = 1;
0256     }
0257 
0258     // pixel quality
0259     const SiPixelQuality* badPixelInfo = &iSetup.getData(siPixelQualityToken_);
0260 
0261     // Tracker global position
0262     const AlignTransform& glbCoord = align::DetectorGlobalPosition(iSetup.getData(gprToken_), DetId(DetId::Tracker));
0263 
0264     // Convert AlignTransform::Translation to GlobalVector using the appropriate constructor
0265     GlobalVector globalTkPosition(glbCoord.translation().x(), glbCoord.translation().y(), glbCoord.translation().z());
0266 
0267     // loop over bclabels
0268     for (const auto& label : bcLabels_) {
0269       // init tree content
0270       PixelBaryCentreAnalyzer::initBC();
0271 
0272       // Get TkAlign from EventSetup:
0273       const Alignments* alignments = &iSetup.getData(tkAlignTokens_[label]);
0274       std::vector<AlignTransform> tkAlignments = alignments->m_align;
0275 
0276       // PIX
0277       GlobalVector barycentre_PIX(0.0, 0.0, 0.0);
0278       // BPIX
0279       GlobalVector barycentre_BPIX(0.0, 0.0, 0.0);
0280       float nmodules_BPIX(0.);
0281       // FPIX
0282       GlobalVector barycentre_FPIX(0.0, 0.0, 0.0);
0283       float nmodules_FPIX(0.);
0284 
0285       // Per-layer/ladder barycentre for BPIX
0286       std::map<int, std::map<int, float>> nmodules_bpix;           // layer-ladder
0287       std::map<int, std::map<int, GlobalVector>> barycentre_bpix;  // layer-ladder
0288 
0289       // Per-disk/ring barycentre for FPIX
0290       std::map<int, std::map<int, float>> nmodules_fpix;           // disk-ring
0291       std::map<int, std::map<int, GlobalVector>> barycentre_fpix;  // disk-ring
0292 
0293       // Loop over tracker module
0294       for (const auto& ali : tkAlignments) {
0295         //DetId
0296         const DetId& detId = DetId(ali.rawId());
0297         // remove bad module
0298         if (usePixelQuality_ && badPixelInfo->IsModuleBad(detId))
0299           continue;
0300 
0301         // alignment for a given module
0302         GlobalVector ali_translation(ali.translation().x(), ali.translation().y(), ali.translation().z());
0303 
0304         int subid = DetId(detId).subdetId();
0305         // BPIX
0306         if (subid == PixelSubdetector::PixelBarrel) {
0307           nmodules_BPIX += 1;
0308           barycentre_BPIX += ali_translation;
0309           barycentre_PIX += ali_translation;
0310 
0311           int layer = tkTopo->pxbLayer(detId);
0312           int ladder = tkTopo->pxbLadder(detId);
0313           nmodules_bpix[layer][ladder] += 1;
0314           barycentre_bpix[layer][ladder] += ali_translation;
0315 
0316         }  // BPIX
0317 
0318         // FPIX
0319         if (subid == PixelSubdetector::PixelEndcap) {
0320           nmodules_FPIX += 1;
0321           barycentre_FPIX += ali_translation;
0322           barycentre_PIX += ali_translation;
0323 
0324           int disk = tkTopo->pxfDisk(detId);
0325           int quadrant = PixelEndcapName(detId, tkTopo, phase_).halfCylinder();
0326           if (quadrant < 3)
0327             disk *= -1;
0328 
0329           int ring = -9999;
0330           if (phase_ == 0) {
0331             ring = 1 + (tkTopo->pxfPanel(detId) + tkTopo->pxfModule(detId.rawId()) > 3);
0332           } else if (phase_ == 1) {
0333             ring = PixelEndcapName(detId, tkTopo, phase_).ringName();
0334           }
0335 
0336           nmodules_fpix[disk][ring] += 1;
0337           barycentre_fpix[disk][ring] += ali_translation;
0338 
0339         }  // FPIX
0340 
0341       }  // loop over tracker module
0342 
0343       //PIX
0344       float nmodules_PIX = nmodules_BPIX + nmodules_FPIX;
0345       barycentre_PIX *= (1.0 / nmodules_PIX);
0346       barycentre_PIX += globalTkPosition;
0347       PIX_ = GlobalPoint(barycentre_PIX.x(), barycentre_PIX.y(), barycentre_PIX.z());
0348       vPIX_ = SimplePoint(PIX_);
0349 
0350       //BPIX
0351       barycentre_BPIX *= (1.0 / nmodules_BPIX);
0352       barycentre_BPIX += globalTkPosition;
0353       BPIX_ = GlobalPoint(barycentre_BPIX.x(), barycentre_BPIX.y(), barycentre_BPIX.z());
0354       vBPIX_ = SimplePoint(BPIX_);
0355       //FPIX
0356       barycentre_FPIX *= (1.0 / nmodules_FPIX);
0357       barycentre_FPIX += globalTkPosition;
0358       FPIX_ = GlobalPoint(barycentre_FPIX.x(), barycentre_FPIX.y(), barycentre_FPIX.z());
0359       vFPIX_ = SimplePoint(FPIX_);
0360       // Pixel substructures
0361 
0362       // BPix barycentre per-layer/per-ladder
0363       // assuming each ladder has the same number of modules in the same layer
0364       // inner =  flipped; outer = non-flipped
0365       //
0366       // Phase 0: Outer ladders are odd for layer 1,3 and even for layer 2
0367       // Phase 1: Outer ladders are odd for layer 4 and even for layer 1,2,3
0368       //
0369 
0370       int nmodules_BPIX_Flipped = 0;
0371       int nmodules_BPIX_NonFlipped = 0;
0372       GlobalVector BPIX_Flipped(0.0, 0.0, 0.0);
0373       GlobalVector BPIX_NonFlipped(0.0, 0.0, 0.0);
0374 
0375       // loop over layers
0376       for (const auto& il : barycentre_bpix) {
0377         int layer = il.first;
0378 
0379         int nmodulesLayer = 0;
0380         int nmodulesLayer_Flipped = 0;
0381         int nmodulesLayer_NonFlipped = 0;
0382         GlobalVector BPIXLayer(0.0, 0.0, 0.0);
0383         GlobalVector BPIXLayer_Flipped(0.0, 0.0, 0.0);
0384         GlobalVector BPIXLayer_NonFlipped(0.0, 0.0, 0.0);
0385 
0386         // loop over ladder
0387         std::map<int, GlobalVector> barycentreLayer = barycentre_bpix[layer];
0388         for (const auto& it : barycentreLayer) {
0389           int ladder = it.first;
0390           //BPIXLayerLadder_[layer][ladder] = (1.0/nmodules[layer][ladder])*barycentreLayer[ladder] + globalTkPosition;
0391 
0392           nmodulesLayer += nmodules_bpix[layer][ladder];
0393           BPIXLayer += barycentreLayer[ladder];
0394 
0395           // Phase-1
0396           //
0397           // Phase 1: Outer ladders are odd for layer 4 and even for layer 1,2,3
0398           if (phase_ == 1) {
0399             if (layer != 4) {  // layer 1-3
0400 
0401               if (ladder % 2 != 0) {  // odd ladder = outer ladder = unflipped
0402                 nmodulesLayer_NonFlipped += nmodules_bpix[layer][ladder];
0403                 BPIXLayer_NonFlipped += barycentreLayer[ladder];
0404               } else {  // even ladder = inner ladder = flipped
0405                 nmodulesLayer_Flipped += nmodules_bpix[layer][ladder];
0406                 BPIXLayer_Flipped += barycentreLayer[ladder];
0407               }
0408             } else {  // layer-4
0409 
0410               if (ladder % 2 != 0) {  // odd ladder = inner = flipped
0411                 nmodulesLayer_Flipped += nmodules_bpix[layer][ladder];
0412                 BPIXLayer_Flipped += barycentreLayer[ladder];
0413               } else {  //even ladder = outer ladder  = unflipped
0414                 nmodulesLayer_NonFlipped += nmodules_bpix[layer][ladder];
0415                 BPIXLayer_NonFlipped += barycentreLayer[ladder];
0416               }
0417             }
0418 
0419           }  // phase-1
0420 
0421           // Phase-0
0422           //
0423           // Phase 0: Outer ladders are odd for layer 1,3 and even for layer 2
0424           if (phase_ == 0) {
0425             if (layer == 2) {  // layer-2
0426 
0427               if (ladder % 2 != 0) {  // odd ladder = inner = flipped
0428                 nmodulesLayer_Flipped += nmodules_bpix[layer][ladder];
0429                 BPIXLayer_Flipped += barycentreLayer[ladder];
0430               } else {
0431                 nmodulesLayer_NonFlipped += nmodules_bpix[layer][ladder];
0432                 BPIXLayer_NonFlipped += barycentreLayer[ladder];
0433               }
0434             } else {  // layer-1,3
0435 
0436               if (ladder % 2 == 0) {  // even ladder = inner = flipped
0437                 nmodulesLayer_Flipped += nmodules_bpix[layer][ladder];
0438                 BPIXLayer_Flipped += barycentreLayer[ladder];
0439               } else {  // odd ladder = outer = non-flipped
0440                 nmodulesLayer_NonFlipped += nmodules_bpix[layer][ladder];
0441                 BPIXLayer_NonFlipped += barycentreLayer[ladder];
0442               }
0443             }
0444 
0445           }  // phase-0
0446 
0447         }  //loop over ladders
0448 
0449         // total BPIX flipped/non-flipped
0450         BPIX_Flipped += BPIXLayer_Flipped;
0451         BPIX_NonFlipped += BPIXLayer_NonFlipped;
0452         nmodules_BPIX_Flipped += nmodulesLayer_Flipped;
0453         nmodules_BPIX_NonFlipped += nmodulesLayer_NonFlipped;
0454 
0455         //BPIX per-layer
0456         BPIXLayer *= (1.0 / nmodulesLayer);
0457         BPIXLayer += globalTkPosition;
0458         BPIXLayer_Flipped *= (1.0 / nmodulesLayer_Flipped);
0459         BPIXLayer_Flipped += globalTkPosition;
0460         BPIXLayer_NonFlipped *= (1.0 / nmodulesLayer_NonFlipped);
0461         BPIXLayer_NonFlipped += globalTkPosition;
0462 
0463         BPIXLayer_[layer - 1] = GlobalPoint(BPIXLayer.x(), BPIXLayer.y(), BPIXLayer.z());
0464         vBPIXLayer_[layer - 1] = SimplePoint(BPIXLayer_[layer - 1]);
0465         BPIXLayer_Flipped_[layer - 1] =
0466             GlobalPoint(BPIXLayer_Flipped.x(), BPIXLayer_Flipped.y(), BPIXLayer_Flipped.z());
0467         vBPIXLayer_Flipped_[layer - 1] = SimplePoint(BPIXLayer_Flipped_[layer - 1]);
0468         BPIXLayer_NonFlipped_[layer - 1] =
0469             GlobalPoint(BPIXLayer_NonFlipped.x(), BPIXLayer_NonFlipped.y(), BPIXLayer_NonFlipped.z());
0470         vBPIXLayer_NonFlipped_[layer - 1] = SimplePoint(BPIXLayer_NonFlipped_[layer - 1]);
0471         BPIXLayer_DiffFlippedNonFlipped_[layer - 1] = GlobalPoint(BPIXLayer_Flipped.x() - BPIXLayer_NonFlipped.x(),
0472                                                                   BPIXLayer_Flipped.y() - BPIXLayer_NonFlipped.y(),
0473                                                                   BPIXLayer_Flipped.z() - BPIXLayer_NonFlipped.z());
0474         vBPIXLayer_DiffFlippedNonFlipped_[layer - 1] = SimplePoint(BPIXLayer_DiffFlippedNonFlipped_[layer - 1]);
0475 
0476       }  // loop over layers
0477 
0478       BPIX_Flipped *= (1.0 / nmodules_BPIX_Flipped);
0479       BPIX_Flipped += globalTkPosition;
0480       BPIX_Flipped_ = GlobalPoint(BPIX_Flipped.x(), BPIX_Flipped.y(), BPIX_Flipped.z());
0481       vBPIX_Flipped_ = SimplePoint(BPIX_Flipped_);
0482       BPIX_NonFlipped *= (1.0 / nmodules_BPIX_NonFlipped);
0483       BPIX_NonFlipped += globalTkPosition;
0484       BPIX_NonFlipped_ = GlobalPoint(BPIX_NonFlipped.x(), BPIX_NonFlipped.y(), BPIX_NonFlipped.z());
0485       vBPIX_NonFlipped_ = SimplePoint(BPIX_NonFlipped_);
0486       BPIX_DiffFlippedNonFlipped_ = GlobalPoint(BPIX_Flipped.x() - BPIX_NonFlipped.x(),
0487                                                 BPIX_Flipped.y() - BPIX_NonFlipped.y(),
0488                                                 BPIX_Flipped.z() - BPIX_NonFlipped.z());
0489       vBPIX_DiffFlippedNonFlipped_ = SimplePoint(BPIX_DiffFlippedNonFlipped_);
0490 
0491       // FPIX substructures per-(signed)disk/per-ring
0492       int nmodules_FPIX_plus = 0;
0493       int nmodules_FPIX_minus = 0;
0494       GlobalVector FPIX_plus(0.0, 0.0, 0.0);
0495       GlobalVector FPIX_minus(0.0, 0.0, 0.0);
0496       // loop over disks
0497 
0498       for (const auto& id : barycentre_fpix) {
0499         int disk = id.first;
0500 
0501         int nmodulesDisk = 0;
0502         GlobalVector FPIXDisk(0.0, 0.0, 0.0);
0503 
0504         std::map<int, GlobalVector> baryCentreDisk = id.second;
0505         for (const auto& ir : baryCentreDisk) {
0506           int ring = ir.first;
0507           nmodulesDisk += nmodules_fpix[disk][ring];
0508           FPIXDisk += ir.second;
0509           if (disk > 0) {
0510             nmodules_FPIX_plus += nmodules_fpix[disk][ring];
0511             FPIX_plus += ir.second;
0512           }
0513           if (disk < 0) {
0514             nmodules_FPIX_minus += nmodules_fpix[disk][ring];
0515             FPIX_minus += ir.second;
0516           }
0517 
0518         }  // loop over rings
0519 
0520         FPIXDisk *= (1.0 / nmodulesDisk);
0521         FPIXDisk += globalTkPosition;
0522 
0523         if (disk > 0) {
0524           FPIXDisks_plus_[disk - 1] = GlobalPoint(FPIXDisk.x(), FPIXDisk.y(), FPIXDisk.z());
0525           vFPIXDisks_plus_[disk - 1] = SimplePoint(FPIXDisks_plus_[disk - 1]);
0526         }
0527         if (disk < 0) {
0528           FPIXDisks_minus_[-disk - 1] = GlobalPoint(FPIXDisk.x(), FPIXDisk.y(), FPIXDisk.z());
0529           vFPIXDisks_minus_[-disk - 1] = SimplePoint(FPIXDisks_minus_[-disk - 1]);
0530         }
0531       }  // loop over disks
0532 
0533       FPIX_plus *= (1.0 / nmodules_FPIX_plus);
0534       FPIX_plus += globalTkPosition;
0535       FPIX_plus_ = GlobalPoint(FPIX_plus.x(), FPIX_plus.y(), FPIX_plus.z());
0536       vFPIX_plus_ = SimplePoint(FPIX_plus_);
0537       FPIX_minus *= (1.0 / nmodules_FPIX_minus);
0538       FPIX_minus += globalTkPosition;
0539       FPIX_minus_ = GlobalPoint(FPIX_minus.x(), FPIX_minus.y(), FPIX_minus.z());
0540       vFPIX_minus_ = SimplePoint(FPIX_minus_);
0541 
0542       bcTrees_[label]->Fill();
0543 
0544     }  // bcLabels_
0545 
0546   }  // check for new IOV for TKAlign
0547 
0548   // beamspot
0549   if (prepareBS) {
0550     // loop over bsLabels_
0551     for (const auto& label : bsLabels_) {
0552       // init bstree content
0553       PixelBaryCentreAnalyzer::initBS();
0554 
0555       // Get BeamSpot from EventSetup
0556       const BeamSpotObjects* mybeamspot = &iSetup.getData(bsTokens_[label]);
0557 
0558       BS_ = GlobalPoint(mybeamspot->x(), mybeamspot->y(), mybeamspot->z());
0559       vBS_ = SimplePoint(BS_);
0560 
0561       bsTrees_[label]->Fill();
0562     }  // bsLabels_
0563 
0564   }  // check for new IOV for BS
0565 }
0566 
0567 // ------------ method called once each job just before starting event loop  ------------
0568 void PixelBaryCentreAnalyzer::beginJob() {
0569   // init bc bs trees
0570   for (const auto& label : bsLabels_) {
0571     std::string treeName = "BeamSpot";
0572     if (!label.empty())
0573       treeName = "BeamSpot_";
0574     treeName += label;
0575 
0576     bsTrees_[label] = tFileService->make<TTree>(TString(treeName), "PixelBarycentre analyzer ntuple");
0577 
0578     bsTrees_[label]->Branch("run", &run_, "run/I");
0579     bsTrees_[label]->Branch("ls", &ls_, "ls/I");
0580 
0581     bsTrees_[label]->Branch("BS", &vBS_, "x/F:y/F:z/F");
0582 
0583   }  // bsLabels_
0584 
0585   for (const auto& label : bcLabels_) {
0586     std::string treeName = "PixelBarycentre";
0587     if (!label.empty())
0588       treeName = "PixelBarycentre_";
0589     treeName += label;
0590     bcTrees_[label] = tFileService->make<TTree>(TString(treeName), "PixelBarycentre analyzer ntuple");
0591 
0592     bcTrees_[label]->Branch("run", &run_, "run/I");
0593     bcTrees_[label]->Branch("ls", &ls_, "ls/I");
0594 
0595     bcTrees_[label]->Branch("PIX", &vPIX_, "x/F:y/F:z/F");
0596 
0597     bcTrees_[label]->Branch("BPIX", &vBPIX_, "x/F:y/F:z/F");
0598     bcTrees_[label]->Branch("BPIX_Flipped", &vBPIX_Flipped_, "x/F:y/F:z/F");
0599     bcTrees_[label]->Branch("BPIX_NonFlipped", &vBPIX_NonFlipped_, "x/F:y/F:z/F");
0600     bcTrees_[label]->Branch("BPIX_DiffFlippedNonFlipped", &vBPIX_DiffFlippedNonFlipped_, "x/F:y/F:z/F");
0601 
0602     bcTrees_[label]->Branch("FPIX", &vFPIX_, "x/F:y/F:z/F");
0603     bcTrees_[label]->Branch("FPIX_plus", &vFPIX_plus_, "x/F:y/F:z/F");
0604     bcTrees_[label]->Branch("FPIX_minus", &vFPIX_minus_, "x/F:y/F:z/F");
0605 
0606     //per-layer
0607     for (unsigned int i = 0; i < nPixelLayers; i++) {
0608       TString structure = "BPIXLYR";
0609       int layer = i + 1;
0610       structure += layer;
0611 
0612       bcTrees_[label]->Branch(structure, &vBPIXLayer_[i], "x/F:y/F:z/F");
0613       bcTrees_[label]->Branch(structure + "_Flipped", &vBPIXLayer_Flipped_[i], "x/F:y/F:z/F");
0614       bcTrees_[label]->Branch(structure + "_NonFlipped", &vBPIXLayer_NonFlipped_[i], "x/F:y/F:z/F");
0615       bcTrees_[label]->Branch(
0616           structure + "_DiffFlippedNonFlipped", &vBPIXLayer_DiffFlippedNonFlipped_[i], "x/F:y/F:z/F");
0617     }
0618 
0619     //per-disk/ring
0620     for (unsigned int i = 0; i < nPixelDiscs; i++) {
0621       TString structure = "FPIXDisk_plus";
0622       int disk = i + 1;
0623       structure += disk;
0624       bcTrees_[label]->Branch(structure, &vFPIXDisks_plus_[i], "x/F:y/F:z/F");
0625 
0626       structure = "FPIXDisk_minus";
0627       structure += disk;
0628       bcTrees_[label]->Branch(structure, &vFPIXDisks_minus_[i], "x/F:y/F:z/F");
0629     }
0630 
0631   }  // bcLabels_
0632 }
0633 
0634 // ------------ method called once each job just after ending the event loop  ------------
0635 void PixelBaryCentreAnalyzer::endJob() {
0636   bcLabels_.clear();
0637   bsLabels_.clear();
0638 
0639   bcTrees_.clear();
0640   bsTrees_.clear();
0641 }
0642 
0643 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0644 void PixelBaryCentreAnalyzer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0645   edm::ParameterSetDescription desc;
0646   desc.setComment("Validates alignment payloads by providing the position of the pixel barycenter positions");
0647   desc.addUntracked<bool>("usePixelQuality", false);
0648   desc.addUntracked<std::vector<std::string>>("tkAlignLabels", {});
0649   desc.addUntracked<std::vector<std::string>>("beamSpotLabels", {});
0650   descriptions.addWithDefaultLabel(desc);
0651 }
0652 
0653 //define this as a plug-in
0654 DEFINE_FWK_MODULE(PixelBaryCentreAnalyzer);