Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:15:30

0001 // -*- C++ -*-
0002 //
0003 /* 
0004  Description: <one line class summary>
0005 
0006  Implementation:
0007      <Notes on implementation>
0008 */
0009 
0010 //
0011 // Original Author:  Riccardo Ranieri
0012 //         Created:  Tue Feb 27 22:22:22 CEST 2007
0013 //
0014 //
0015 
0016 // system include files
0017 #include <memory>
0018 
0019 // user include files
0020 #include "FWCore/Framework/interface/Frameworkfwd.h"
0021 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0022 
0023 #include "FWCore/Framework/interface/Event.h"
0024 #include "FWCore/Framework/interface/EventSetup.h"
0025 #include "FWCore/Framework/interface/ESHandle.h"
0026 #include "FWCore/Framework/interface/MakerMacros.h"
0027 
0028 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0029 #include "Geometry/CommonDetUnit/interface/TrackingGeometry.h"
0030 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0031 #include "Geometry/TrackerNumberingBuilder/interface/GeometricDet.h"
0032 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0033 #include "Geometry/CommonTopologies/interface/StripTopology.h"
0034 #include "Geometry/CommonDetUnit/interface/PixelGeomDetType.h"
0035 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetType.h"
0036 
0037 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0038 #include "DataFormats/GeometrySurface/interface/BoundSurface.h"
0039 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0040 
0041 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0042 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0043 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0044 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0045 #include "Geometry/TrackerNumberingBuilder/interface/CmsTrackerStringToEnum.h"
0046 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0047 
0048 // output
0049 #include <iostream>
0050 #include <fstream>
0051 #include <iomanip>
0052 #include <cmath>
0053 #include <bitset>
0054 //
0055 
0056 //
0057 // class decleration
0058 //
0059 
0060 //double PI = 3.141592654;
0061 
0062 class ModuleNumbering : public edm::one::EDAnalyzer<> {
0063 public:
0064   explicit ModuleNumbering(const edm::ParameterSet&);
0065   ~ModuleNumbering() override;
0066 
0067   void beginJob() override {}
0068   void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0069   void endJob() override {}
0070 
0071 private:
0072   // ----------member data ---------------------------
0073   void fillModuleVariables(const GeometricDet* module, double& polarRadius, double& phiRad, double& z);
0074   double changePhiRange_From_ZeroTwoPi_To_MinusPiPlusPi(double phiRad);
0075   double changePhiRange_From_MinusPiPlusPi_To_MinusTwoPiZero(double phiRad);
0076   double changePhiRange_From_MinusPiPlusPi_To_ZeroTwoPi(double phiRad);
0077 
0078   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> tokTopo_;
0079   const edm::ESGetToken<GeometricDet, IdealGeometryRecord> tokGeo_;
0080 
0081   //
0082   // counters
0083   unsigned int iOK;
0084   unsigned int iERROR;
0085   //
0086 };
0087 
0088 //
0089 // constants, enums and typedefs
0090 //
0091 
0092 //
0093 // static data member definitions
0094 //
0095 static const double tolerance_space = 1.000;  // 1.000 mm
0096 static const double tolerance_angle = 0.001;  // 0.001 rad
0097 
0098 //
0099 // constructors and destructor
0100 //
0101 ModuleNumbering::ModuleNumbering(const edm::ParameterSet& iConfig)
0102     : tokTopo_(esConsumes<TrackerTopology, TrackerTopologyRcd>()),
0103       tokGeo_(esConsumes<GeometricDet, IdealGeometryRecord>()) {
0104   //now do what ever initialization is needed
0105 }
0106 
0107 ModuleNumbering::~ModuleNumbering() {
0108   // do anything here that needs to be done at desctruction time
0109   // (e.g. close files, deallocate resources etc.)
0110 }
0111 
0112 //
0113 // member functions
0114 //
0115 
0116 void ModuleNumbering::fillModuleVariables(const GeometricDet* module, double& polarRadius, double& phiRad, double& z) {
0117   // module variables
0118   polarRadius = std::sqrt(module->translation().X() * module->translation().X() +
0119                           module->translation().Y() * module->translation().Y());
0120   phiRad = atan2(module->translation().Y(), module->translation().X());
0121   // tolerance near phi=0
0122   if (std::abs(phiRad) < tolerance_angle)
0123     phiRad = 0.0;
0124   // negative phi: from [-PI,+PI) to [0,2PI)
0125   if (phiRad < 0)
0126     phiRad += 2 * M_PI;
0127   //
0128   z = module->translation().Z();
0129   //
0130 }
0131 
0132 double ModuleNumbering::changePhiRange_From_ZeroTwoPi_To_MinusPiPlusPi(double phiRad) {
0133   double new_phiRad = phiRad;
0134   // tolerance near phi=PI
0135   if (std::abs(new_phiRad - M_PI) < tolerance_angle)
0136     new_phiRad = M_PI;
0137   // phi greater than PI: from [0,2PI) to [-PI,+PI)
0138   if (new_phiRad > M_PI)
0139     new_phiRad -= 2 * M_PI;
0140   //
0141   return new_phiRad;
0142 }
0143 
0144 double ModuleNumbering::changePhiRange_From_MinusPiPlusPi_To_MinusTwoPiZero(double phiRad) {
0145   double new_phiRad = phiRad;
0146   // tolerance near phi=PI
0147   if (std::abs(std::abs(new_phiRad) - M_PI) < tolerance_angle)
0148     new_phiRad = M_PI;
0149   // phi greater than PI: from [-PI,+PI) to [0,2PI)
0150   if (new_phiRad > 0)
0151     new_phiRad -= 2 * M_PI;
0152   //
0153   return new_phiRad;
0154 }
0155 
0156 double ModuleNumbering::changePhiRange_From_MinusPiPlusPi_To_ZeroTwoPi(double phiRad) {
0157   double new_phiRad = phiRad;
0158   // tolerance near phi=PI
0159   if (std::abs(std::abs(new_phiRad) - M_PI) < tolerance_angle)
0160     new_phiRad = M_PI;
0161   // phi greater than PI: from [-PI,+PI) to [0,2PI)
0162   if (new_phiRad < 0)
0163     new_phiRad += 2 * M_PI;
0164   //
0165   return new_phiRad;
0166 }
0167 
0168 // ------------ method called to produce the data  ------------
0169 void ModuleNumbering::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0170   //Retrieve tracker topology from geometry
0171   const TrackerTopology* tTopo = &iSetup.getData(tokTopo_);
0172 
0173   edm::LogInfo("ModuleNumbering") << "begins";
0174 
0175   // reset counters
0176   iOK = 0;
0177   iERROR = 0;
0178   //
0179 
0180   //
0181   // get the GeometricDet
0182   //
0183   auto const& rDD = iSetup.getHandle(tokGeo_);
0184   edm::LogInfo("ModuleNumbering") << " Top node is  " << rDD.product() << " " << rDD.product()->name() << std::endl;
0185   edm::LogInfo("ModuleNumbering") << "    radLength " << rDD.product()->radLength() << "\n"
0186                                   << "           xi " << rDD.product()->xi() << "\n"
0187                                   << " PixelROCRows " << rDD.product()->pixROCRows() << "\n"
0188                                   << "   PixROCCols " << rDD.product()->pixROCCols() << "\n"
0189                                   << "   PixelROC_X " << rDD.product()->pixROCx() << "\n"
0190                                   << "   PixelROC_Y " << rDD.product()->pixROCy() << "\n"
0191                                   << "TrackerStereoDetectors " << (rDD.product()->stereo() ? "true" : "false") << "\n"
0192                                   << "SiliconAPVNumber " << rDD.product()->siliconAPVNum() << "\n";
0193   std::vector<int> nv = rDD.product()->navType();
0194   edm::LogInfo("ModuleNumbering").log([&](auto& log) {
0195     for (auto it : nv)
0196       log << it << ", ";
0197   });
0198 
0199   edm::LogInfo("ModuleNumbering") << " And Contains  Daughters: " << rDD.product()->deepComponents().size()
0200                                   << std::endl;
0201   // output file
0202   const std::string& outputFileName =
0203       (!rDD.product()->isFromDD4hep() ? "ModuleNumbering.log" : "ModuleNumbering_dd4hep.log");
0204   std::ofstream Output(outputFileName, std::ios::out);
0205 
0206   //
0207   //first instance tracking geometry
0208   //
0209 
0210   std::vector<const GeometricDet*> modules = (*rDD).deepComponents();
0211   std::map<uint32_t, const GeometricDet*> mapDetIdToGeometricDet;
0212 
0213   for (auto& module : modules) {
0214     mapDetIdToGeometricDet[module->geographicalId().rawId()] = module;
0215   }
0216 
0217   // Debug variables
0218   //
0219   uint32_t myDetId = 0;
0220   unsigned int iDetector = 1;
0221   unsigned int nSubDetectors = 6;
0222   //
0223   double polarRadius = 0.0;
0224   double phiRad = 0.0;
0225   double z = 0.0;
0226   //
0227 
0228   Output << "************************ List of modules with positions ************************" << std::endl;
0229   Output << std::fixed << std::setprecision(4);  // set as default 4 decimal digits (0.1 um or 0.1 rad)
0230 
0231   for (unsigned int iSubDetector = 1; iSubDetector <= nSubDetectors; iSubDetector++) {
0232     // modules
0233     switch (iSubDetector) {
0234         // PXB
0235       case 1: {
0236         break;
0237       }
0238 
0239         // PXF
0240       case 2: {
0241         break;
0242       }
0243 
0244         // TIB
0245       case 3: {
0246         // TIB loop
0247         // number of strings per layer 1/4 int/ext
0248         unsigned int string_int_ext_TIB[8] = {26, 30, 34, 38, 44, 46, 52, 56};
0249         unsigned int mod_type_TIB[8] = {1, 2, 1, 2, 0, 0, 0, 0};  // first and last type for module type loop
0250         unsigned int nLayers = 4;
0251         unsigned int nModules = 3;
0252         // debug variables
0253         double layer_R = 0.0;
0254         double layer_R_previous = 0.0;
0255         double side_z = 0.0;
0256         double side_z_previous = -10000.0;
0257         double part_R = 0.0;
0258         double part_R_previous = 0.0;
0259         double string_phi = 0.0;
0260         double string_phi_previous = 0.0;
0261         double module_z = 0.0;
0262         double module_z_previous = 0.0;
0263         //
0264         for (unsigned int iLayer = 1; iLayer <= nLayers; iLayer++) {  // Layer: 1,...,nLayers
0265           for (unsigned int iSide = 1; iSide <= 2; iSide++) {         // Side: 1:- 2:+
0266             for (unsigned int iPart = 1; iPart <= 2; iPart++) {       // Part: 1:int 2:ext
0267               unsigned int jString = (2 * (iLayer - 1)) + (iPart - 1);
0268               for (unsigned int iString = 1; iString <= string_int_ext_TIB[jString];
0269                    iString++) {                                                   // String: 1,...,nStrings
0270                 for (unsigned int iModule = 1; iModule <= nModules; iModule++) {  // Module: 1,...,nModules
0271                   for (unsigned int iType = mod_type_TIB[2 * (iLayer - 1)]; iType <= mod_type_TIB[2 * (iLayer - 1) + 1];
0272                        iType++) {  // Module Type: 0 (ss) 1-2 (ds stereo and rphi)
0273 
0274                     myDetId = 0;
0275                     // detector
0276                     myDetId <<= 4;
0277                     myDetId |= iDetector;
0278                     // subdetector
0279                     myDetId <<= 3;
0280                     myDetId |= iSubDetector;
0281                     // not used
0282                     myDetId <<= 8;
0283                     myDetId |= 0;
0284                     // layer
0285                     myDetId <<= 3;
0286                     myDetId |= iLayer;
0287                     // side
0288                     myDetId <<= 2;
0289                     myDetId |= iSide;
0290                     // part
0291                     myDetId <<= 2;
0292                     myDetId |= iPart;
0293                     // string number
0294                     myDetId <<= 6;
0295                     myDetId |= iString;
0296                     // module number
0297                     myDetId <<= 2;
0298                     myDetId |= iModule;
0299                     // module type
0300                     myDetId <<= 2;
0301                     myDetId |= iType;
0302                     //
0303                     std::bitset<32> binary_myDetId(myDetId);
0304                     Output << std::endl << std::endl;
0305                     Output << " ******** myDet Id = " << myDetId << " (" << binary_myDetId << ")" << std::endl;
0306                     //
0307                     unsigned int rawid = mapDetIdToGeometricDet[myDetId]->geographicalId().rawId();
0308                     std::bitset<32> binary_detid(rawid);
0309                     GeometricDet::nav_type detNavType = mapDetIdToGeometricDet[myDetId]->navType();
0310                     //
0311                     Output << "            raw Id = " << rawid << " (" << binary_detid << ")"
0312                            << "\t nav type = " << GeometricDet::printNavType(&detNavType.front(), detNavType.size())
0313                            << std::endl;
0314 
0315                     // variables
0316                     fillModuleVariables(mapDetIdToGeometricDet[myDetId], polarRadius, phiRad, z);
0317                     layer_R = polarRadius;
0318                     side_z = z;
0319                     part_R = polarRadius;
0320                     string_phi = phiRad;
0321                     module_z = z;
0322                     //
0323 
0324                     Output << "\t R = " << polarRadius << "\t phi = " << phiRad << "\t z = " << z << std::endl;
0325 
0326                     // Module Info
0327 
0328                     std::string name = mapDetIdToGeometricDet[myDetId]->name();
0329                     unsigned int theLayer = tTopo->tibLayer(rawid);
0330                     std::vector<unsigned int> theString = tTopo->tibStringInfo(rawid);
0331                     unsigned int theModule = tTopo->tibModule(rawid);
0332                     std::string side;
0333                     std::string part;
0334                     side = (theString[0] == 1) ? "-" : "+";
0335                     part = (theString[1] == 1) ? "int" : "ext";
0336                     Output << " TIB" << side << "\t"
0337                            << "Layer " << theLayer << " " << part << "\t"
0338                            << "string " << theString[2] << "\t"
0339                            << " module " << theModule << " " << name << std::endl;
0340                     //
0341 
0342                     // module: |z| check
0343                     Output << "\t # ";
0344                     if ((std::abs(module_z) - std::abs(module_z_previous)) < (0 + tolerance_space)) {
0345                       Output << "\t ERROR |z| ordering not respected in module ";
0346                       iERROR++;
0347                     } else {
0348                       Output << "\t OK"
0349                              << " |z| ordering in module ";
0350                       iOK++;
0351                     }
0352                     Output << iModule - 1 << " to " << iModule << " (" << module_z_previous << " --> " << module_z
0353                            << ")" << std::endl;
0354                     //
0355                   }  // type loop
0356 
0357                   //
0358                   module_z_previous = module_z;
0359                   //
0360                 }  // module loop
0361 
0362                 // string: phi check
0363                 Output << "\t ## ";
0364                 if ((string_phi - string_phi_previous) < (0 - tolerance_angle)) {
0365                   Output << "\t ERROR phi ordering not respected in string ";
0366                   iERROR++;
0367                 } else {
0368                   Output << "\t OK"
0369                          << " phi ordering in string ";
0370                   iOK++;
0371                 }
0372                 Output << iString - 1 << " to " << iString << " (" << string_phi_previous << " --> " << string_phi
0373                        << ")" << std::endl;
0374                 //
0375                 string_phi_previous = string_phi;
0376                 module_z_previous = 0.0;
0377                 //
0378               }  // string loop
0379 
0380               // part: R check
0381               Output << "\t ### ";
0382               if ((part_R - part_R_previous) < (0 + tolerance_space)) {
0383                 Output << "\t ERROR R ordering (part int/ext) not respected in layer ";
0384                 iERROR++;
0385               } else {
0386                 Output << "\t OK"
0387                        << " R ordering (part int/ext) in layer ";
0388                 iOK++;
0389               }
0390               Output << iLayer << " part " << iPart - 1 << " to " << iPart << " (" << part_R_previous << " --> "
0391                      << part_R << ")" << std::endl;
0392               //
0393               part_R_previous = part_R;
0394               string_phi_previous = 0.0;
0395               module_z_previous = 0.0;
0396               //
0397             }  // part loop
0398 
0399             // side: z check
0400             Output << "\t #### ";
0401             if ((side_z - side_z_previous) < (0 + tolerance_space)) {
0402               Output << "\t ERROR z ordering (side -/+) not respected in layer ";
0403               iERROR++;
0404             } else {
0405               Output << "\t OK"
0406                      << " z ordering (side -/+) in layer ";
0407               iOK++;
0408             }
0409             Output << iLayer << " side " << iSide - 1 << " to " << iSide << " (" << side_z_previous << " --> " << side_z
0410                    << ")" << std::endl;
0411             //
0412             side_z_previous = side_z;
0413             part_R_previous = 0.0;
0414             string_phi_previous = 0.0;
0415             module_z_previous = 0.0;
0416             //
0417           }  // side loop
0418 
0419           // layer: R check
0420           Output << "\t ##### ";
0421           if ((layer_R - layer_R_previous) < (0 + tolerance_space)) {
0422             Output << "\t ERROR R ordering not respected from layer ";
0423             iERROR++;
0424           } else {
0425             Output << "\t OK"
0426                    << " R ordering in layer ";
0427             iOK++;
0428           }
0429           Output << iLayer - 1 << " to " << iLayer << " (" << layer_R_previous << " --> " << layer_R << ")"
0430                  << std::endl;
0431           //
0432           layer_R_previous = layer_R;
0433           side_z_previous = -10000.0;
0434           part_R_previous = 0.0;
0435           string_phi_previous = 0.0;
0436           module_z_previous = 0.0;
0437           //
0438         }  // layer loop
0439 
0440         break;
0441       }
0442 
0443         // TID
0444       case 4: {
0445         // TID loop
0446         unsigned int modules_TID[3] = {12, 12, 20};               // number of modules per disk
0447         unsigned int mod_type_TID[8] = {1, 2, 1, 2, 0, 0, 0, 0};  // first and last type for module type loop
0448         unsigned int nDisks = 3;
0449         unsigned int nRings = 3;
0450         // debug variables
0451         double side_z = 0.0;
0452         double side_z_previous = 10000.0;
0453         double disk_z = 0.0;
0454         double disk_z_previous = 0.0;
0455         double ring_R = 0.0;
0456         double ring_R_previous = 0.0;
0457         double part_z = 0.0;
0458         double part_z_previous = 0.0;
0459         double module_phi = 0.0;
0460         double module_phi_previous = -M_PI;
0461         //
0462         for (unsigned int iSide = 2; iSide >= 1; iSide--) {           // Side: 1:- 2:+
0463           for (unsigned int iDisk = 1; iDisk <= nDisks; iDisk++) {    // Disk: 1,...,nDisks
0464             for (unsigned int iRing = 1; iRing <= nRings; iRing++) {  // Ring: 1,...,nRings
0465               for (int iPart = 2; iPart >= 1; iPart--) {              // Part: 1:back 2:front
0466                 for (unsigned int iModule = 1; iModule <= modules_TID[iRing - 1];
0467                      iModule++) {  // Module: 1,...,modules in ring
0468                   for (unsigned int iType = mod_type_TID[2 * (iRing - 1)]; iType <= mod_type_TID[2 * (iRing - 1) + 1];
0469                        iType++) {  // Module Type: 0 (ss) 1-2 (ds stereo and rphi)
0470 
0471                     myDetId = 0;
0472                     // detector
0473                     myDetId <<= 4;
0474                     myDetId |= iDetector;
0475                     // subdetector
0476                     myDetId <<= 3;
0477                     myDetId |= iSubDetector;
0478                     // not used
0479                     myDetId <<= 10;
0480                     myDetId |= 0;
0481                     // side
0482                     myDetId <<= 2;
0483                     myDetId |= iSide;
0484                     // disk number
0485                     myDetId <<= 2;
0486                     myDetId |= iDisk;
0487                     // ring number
0488                     myDetId <<= 2;
0489                     myDetId |= iRing;
0490                     // part
0491                     myDetId <<= 2;
0492                     myDetId |= iPart;
0493                     // module number
0494                     myDetId <<= 5;
0495                     myDetId |= iModule;
0496                     // module type
0497                     myDetId <<= 2;
0498                     myDetId |= iType;
0499                     //
0500                     std::bitset<32> binary_myDetId(myDetId);
0501                     Output << std::endl << std::endl;
0502                     Output << " ******** myDet Id = " << myDetId << " (" << binary_myDetId << ")" << std::endl;
0503                     //
0504                     unsigned int rawid = mapDetIdToGeometricDet[myDetId]->geographicalId().rawId();
0505                     std::bitset<32> binary_detid(rawid);
0506                     GeometricDet::nav_type detNavType = mapDetIdToGeometricDet[myDetId]->navType();
0507                     //
0508                     Output << "            raw Id = " << rawid << " (" << binary_detid << ")"
0509                            << "\t nav type = " << GeometricDet::printNavType(&detNavType.front(), detNavType.size())
0510                            << std::endl;
0511 
0512                     // variables
0513                     fillModuleVariables(mapDetIdToGeometricDet[myDetId], polarRadius, phiRad, z);
0514                     side_z = z;
0515                     disk_z = z;
0516                     ring_R = polarRadius;
0517                     part_z = z;
0518                     module_phi = phiRad;
0519                     //
0520 
0521                     Output << "\t R = " << polarRadius << "\t phi = " << phiRad << "\t z = " << z << std::endl;
0522 
0523                     // Module Info
0524 
0525                     std::string name = mapDetIdToGeometricDet[myDetId]->name();
0526                     unsigned int theDisk = tTopo->tidWheel(rawid);
0527                     unsigned int theRing = tTopo->tidRing(rawid);
0528                     std::vector<unsigned int> theModule = tTopo->tidModuleInfo(rawid);
0529                     std::string side;
0530                     std::string part;
0531                     side = (tTopo->tidSide(rawid) == 1) ? "-" : "+";
0532                     part = (theModule[0] == 1) ? "back" : "front";
0533                     Output << " TID" << side << "\t"
0534                            << "Disk " << theDisk << " Ring " << theRing << " " << part << "\t"
0535                            << " module " << theModule[1] << "\t" << name << std::endl;
0536                     //
0537 
0538                     // module: phi check
0539                     Output << "\t # ";
0540                     if ((module_phi - module_phi_previous) < (0 - tolerance_angle)) {
0541                       Output << "\t ERROR phi ordering not respected in module ";
0542                       iERROR++;
0543                     } else {
0544                       Output << "\t OK"
0545                              << " phi ordering in module ";
0546                       iOK++;
0547                     }
0548                     Output << iModule - 1 << " to " << iModule << " (" << module_phi_previous << " --> " << module_phi
0549                            << ")" << std::endl;
0550                     //
0551                   }  // type loop
0552 
0553                   //
0554                   module_phi_previous = module_phi;
0555                   //
0556                 }  // module loop
0557 
0558                 // part: |z| check
0559                 Output << "\t ## ";
0560                 if ((std::abs(part_z) - std::abs(part_z_previous)) < (0 + tolerance_space)) {
0561                   Output << "\t ERROR |z| ordering (front/back) not respected in ring ";
0562                   iERROR++;
0563                 } else {
0564                   Output << "\t OK"
0565                          << " |z| ordering (front/back) in ring ";
0566                   iOK++;
0567                 }
0568                 Output << iRing << " part " << iPart + 1 << " to " << iPart << " (" << part_z_previous << " --> "
0569                        << part_z << ")" << std::endl;
0570                 //
0571                 part_z_previous = part_z;
0572                 module_phi_previous = -M_PI;
0573                 //
0574               }  // part loop
0575 
0576               // ring: R check
0577               Output << "\t ### ";
0578               if ((ring_R - ring_R_previous) < (0 + tolerance_space)) {
0579                 Output << "\t ERROR R ordering not respected in disk ";
0580                 iERROR++;
0581               } else {
0582                 Output << "\t OK"
0583                        << " R ordering in disk ";
0584                 iOK++;
0585               }
0586               Output << iDisk << " ring " << iRing - 1 << " to " << iRing << " (" << ring_R_previous << " --> "
0587                      << ring_R << ")" << std::endl;
0588               //
0589               ring_R_previous = ring_R;
0590               part_z_previous = 0.0;
0591               module_phi_previous = -M_PI;
0592               //
0593             }  // ring loop
0594 
0595             // disk: |z| check
0596             Output << "\t #### ";
0597             if ((std::abs(disk_z) - std::abs(disk_z_previous)) < (0 + tolerance_space)) {
0598               Output << "\t ERROR |z| ordering not respected in disk ";
0599               iERROR++;
0600             } else {
0601               Output << "\t OK"
0602                      << " |z| ordering in disk ";
0603               iOK++;
0604             }
0605             Output << iDisk - 1 << " to " << iDisk << " (" << disk_z_previous << " --> " << disk_z << ")" << std::endl;
0606             //
0607             disk_z_previous = disk_z;
0608             ring_R_previous = 0.0;
0609             part_z_previous = 0.0;
0610             module_phi_previous = -M_PI;
0611             //
0612           }  // disk loop
0613 
0614           // side: z check
0615           Output << "\t ##### ";
0616           if ((side_z - side_z_previous) > (0 + tolerance_space)) {
0617             Output << "\t ERROR z ordering (side -/+) not respected in TID side ";
0618             iERROR++;
0619           } else {
0620             Output << "\t OK"
0621                    << " z ordering (side -/+) in TID side ";
0622             iOK++;
0623           }
0624           Output << iSide + 1 << " to " << iSide << " (" << side_z_previous << " --> " << side_z << ")" << std::endl;
0625           //
0626           side_z_previous = side_z;
0627           disk_z_previous = 0.0;
0628           ring_R_previous = 0.0;
0629           part_z_previous = 0.0;
0630           module_phi_previous = -M_PI;
0631           //
0632         }  // side loop
0633 
0634         break;
0635       }
0636 
0637         // TOB
0638       case 5: {
0639         // TOB loop
0640         unsigned int rod_TOB[8] = {42, 48, 54, 60, 66, 74};  // number of rods per layer 1/6
0641         unsigned int mod_type_TOB[12] = {
0642             1, 2, 1, 2, 0, 0, 0, 0, 0, 0, 0, 0};  // first and last type for module type loop
0643         unsigned int nLayers = 6;
0644         unsigned int nModules = 6;
0645         // debug variables
0646         double layer_R = 0.0;
0647         double layer_R_previous = 0.0;
0648         double side_z = 0.0;
0649         double side_z_previous = -10000.0;
0650         double rod_phi = 0.0;
0651         double rod_phi_previous = 0.0;
0652         double module_z = 0.0;
0653         double module_z_previous = 0.0;
0654         //
0655         for (unsigned int iLayer = 1; iLayer <= nLayers; iLayer++) {            // Layer: 1,...,nLayers
0656           for (unsigned int iSide = 1; iSide <= 2; iSide++) {                   // Side: 1:- 2:+
0657             for (unsigned int iRod = 1; iRod <= rod_TOB[iLayer - 1]; iRod++) {  // Rod: 1,...,nRods
0658               for (unsigned int iModule = 1; iModule <= nModules; iModule++) {  // Module: 1,...,nModules
0659                 for (unsigned int iType = mod_type_TOB[2 * (iLayer - 1)]; iType <= mod_type_TOB[2 * (iLayer - 1) + 1];
0660                      iType++) {  // Module Type: 0 (ss) 1-2 (ds stereo and rphi)
0661 
0662                   myDetId = 0;
0663                   // detector
0664                   myDetId <<= 4;
0665                   myDetId |= iDetector;
0666                   // subdetector
0667                   myDetId <<= 3;
0668                   myDetId |= iSubDetector;
0669                   // not used
0670                   myDetId <<= 8;
0671                   myDetId |= 0;
0672                   // layer
0673                   myDetId <<= 3;
0674                   myDetId |= iLayer;
0675                   // side
0676                   myDetId <<= 2;
0677                   myDetId |= iSide;
0678                   // rod number
0679                   myDetId <<= 7;
0680                   myDetId |= iRod;
0681                   // module number
0682                   myDetId <<= 3;
0683                   myDetId |= iModule;
0684                   // module type
0685                   myDetId <<= 2;
0686                   myDetId |= iType;
0687                   //
0688                   std::bitset<32> binary_myDetId(myDetId);
0689                   Output << std::endl << std::endl;
0690                   Output << " ******** myDet Id = " << myDetId << " (" << binary_myDetId << ")" << std::endl;
0691                   //
0692                   unsigned int rawid = mapDetIdToGeometricDet[myDetId]->geographicalId().rawId();
0693                   std::bitset<32> binary_detid(rawid);
0694                   GeometricDet::nav_type detNavType = mapDetIdToGeometricDet[myDetId]->navType();
0695                   //
0696                   Output << "            raw Id = " << rawid << " (" << binary_detid << ")"
0697                          << "\t nav type = " << GeometricDet::printNavType(&detNavType.front(), detNavType.size())
0698                          << std::endl;
0699 
0700                   // variables
0701                   fillModuleVariables(mapDetIdToGeometricDet[myDetId], polarRadius, phiRad, z);
0702                   layer_R = polarRadius;
0703                   side_z = z;
0704                   rod_phi = phiRad;
0705                   module_z = z;
0706                   //
0707 
0708                   Output << "\t R = " << polarRadius << "\t phi = " << phiRad << "\t z = " << z << std::endl;
0709 
0710                   // Module Info
0711 
0712                   std::string name = mapDetIdToGeometricDet[myDetId]->name();
0713                   unsigned int theLayer = tTopo->tobLayer(rawid);
0714                   std::vector<unsigned int> theRod = tTopo->tobRodInfo(rawid);
0715                   unsigned int theModule = tTopo->tobModule(rawid);
0716                   std::string side;
0717                   std::string part;
0718                   side = (theRod[0] == 1) ? "-" : "+";
0719                   Output << " TOB" << side << "\t"
0720                          << "Layer " << theLayer << "\t"
0721                          << "rod " << theRod[1] << " module " << theModule << "\t" << name << std::endl;
0722                   //
0723 
0724                   // module: |z| check
0725                   Output << "\t # ";
0726                   if ((std::abs(module_z) - std::abs(module_z_previous)) < (0 + tolerance_space)) {
0727                     Output << "\t ERROR |z| ordering not respected in module ";
0728                     iERROR++;
0729                   } else {
0730                     Output << "\t OK"
0731                            << " |z| ordering in module ";
0732                     iOK++;
0733                   }
0734                   Output << iModule - 1 << " to " << iModule << " (" << module_z_previous << " --> " << module_z << ")"
0735                          << std::endl;
0736                   //
0737                 }  // type loop
0738 
0739                 //
0740                 module_z_previous = module_z;
0741                 //
0742               }  // module loop
0743 
0744               // rod: phi check
0745               Output << "\t ## ";
0746               if ((rod_phi - rod_phi_previous) < (0 - tolerance_angle)) {
0747                 Output << "\t ERROR phi ordering not respected in rod ";
0748                 iERROR++;
0749               } else {
0750                 Output << "\t OK"
0751                        << " phi ordering in rod ";
0752                 iOK++;
0753               }
0754               Output << iRod - 1 << " to " << iRod << " (" << rod_phi_previous << " --> " << rod_phi << ")"
0755                      << std::endl;
0756               //
0757               rod_phi_previous = rod_phi;
0758               module_z_previous = 0.0;
0759               //
0760             }  // rod loop
0761 
0762             // side: z check
0763             Output << "\t ### ";
0764             if ((side_z - side_z_previous) < (0 + tolerance_space)) {
0765               Output << "\t ERROR z ordering (side -/+) not respected in layer ";
0766               iERROR++;
0767             } else {
0768               Output << "\t OK"
0769                      << " z ordering (side -/+) in layer ";
0770               iOK++;
0771             }
0772             Output << iLayer << " side " << iSide - 1 << " to " << iSide << " (" << side_z_previous << " --> " << side_z
0773                    << ")" << std::endl;
0774             //
0775             side_z_previous = side_z;
0776             rod_phi_previous = 0.0;
0777             module_z_previous = 0.0;
0778             //
0779           }  // side loop
0780 
0781           // layer: R check
0782           Output << "\t #### ";
0783           if ((layer_R - layer_R_previous) < (0 + tolerance_space)) {
0784             Output << "\t ERROR R ordering not respected from layer ";
0785             iERROR++;
0786           } else {
0787             Output << "\t OK"
0788                    << " R ordering in layer ";
0789             iOK++;
0790           }
0791           Output << iLayer - 1 << " to " << iLayer << " (" << layer_R_previous << " --> " << layer_R << ")"
0792                  << std::endl;
0793           //
0794           layer_R_previous = layer_R;
0795           side_z_previous = -10000.0;
0796           rod_phi_previous = 0.0;
0797           module_z_previous = 0.0;
0798           //
0799         }  // layer loop
0800 
0801         break;
0802       }
0803 
0804         // TEC
0805       case 6: {
0806         // TEC loop
0807         unsigned int nWheels = 9;
0808         unsigned int nPetals = 8;
0809         unsigned int nRings = 7;
0810         unsigned int first_ring_TEC[9] = {1, 1, 1, 2, 2, 2, 3, 3, 4};  // first ring of the wheel
0811         unsigned int modules_ring_TEC[14] = {
0812             1, 2, 1, 2, 2, 3, 3, 4, 3, 2, 3, 4, 5, 5};  // number of modules ring=1,...,nRings back/front
0813         unsigned int mod_type_TEC[14] = {
0814             1, 2, 1, 2, 0, 0, 0, 0, 1, 2, 0, 0, 0, 0};  // first and last type for module type loop (per ring)
0815         // debug variables
0816         double side_z = 0.0;
0817         double side_z_previous = 10000.0;
0818         double wheel_z = 0.0;
0819         double wheel_z_previous = 0.0;
0820         double part_z = 0.0;
0821         double part_z_previous = 0.0;
0822         double petal_phi = 0.0;
0823         double petal_phi_previous = 0.0;
0824         double ring_R = 0.0;
0825         double ring_R_previous = 0.0;
0826         double module_phi = 0.0;
0827         double module_phi_previous = -M_PI;
0828         //
0829         for (unsigned int iSide = 2; iSide >= 1; iSide--) {  // Side: 1:- 2:+
0830           switch (iSide) {
0831               // TEC+
0832             case 2: {
0833               side_z = 0.0;
0834               side_z_previous = 10000.0;
0835               wheel_z = 0.0;
0836               wheel_z_previous = 0.0;
0837               part_z = 0.0;
0838               part_z_previous = 0.0;
0839               petal_phi = 0.0;
0840               petal_phi_previous = 0.0;
0841               ring_R = 0.0;
0842               ring_R_previous = 0.0;
0843               module_phi = 0.0;
0844               module_phi_previous = -M_PI;
0845               break;
0846             }
0847               // TEC-
0848             case 1: {
0849               wheel_z = 0.0;
0850               wheel_z_previous = 0.0;
0851               part_z = 0.0;
0852               part_z_previous = 0.0;
0853               petal_phi = 0.0;
0854               petal_phi_previous = 0.0;
0855               ring_R = 0.0;
0856               ring_R_previous = 0.0;
0857               module_phi = 0.0;
0858               module_phi_previous = 2 * M_PI;
0859               break;
0860             }
0861             default: {
0862               // do nothing
0863             }
0864           }
0865           //
0866           for (unsigned int iWheel = 1; iWheel <= nWheels; iWheel++) {      // Wheel: 1,...,nWheels
0867             for (int iPart = 2; iPart >= 1; iPart--) {                      // Part: 1:back 2:front
0868               for (unsigned int iPetal = 1; iPetal <= nPetals; iPetal++) {  // Petal: 1,...,nPetals
0869                 for (unsigned int iRing = first_ring_TEC[iWheel - 1]; iRing <= nRings;
0870                      iRing++) {  // Ring: first,...,nRings
0871                   unsigned int nModules = modules_ring_TEC[2 * (iRing - 1) + (iPart - 1)];
0872                   for (unsigned int iModule = 1; iModule <= nModules;
0873                        iModule++) {  // Module: 1,...,modules in ring of petal
0874                     for (unsigned int iType = mod_type_TEC[2 * (iRing - 1)]; iType <= mod_type_TEC[2 * (iRing - 1) + 1];
0875                          iType++) {  // Module Type: 0 (ss) 1-2 (ds stereo and rphi)
0876 
0877                       myDetId = 0;
0878                       // detector
0879                       myDetId <<= 4;
0880                       myDetId |= iDetector;
0881                       // subdetector
0882                       myDetId <<= 3;
0883                       myDetId |= iSubDetector;
0884                       // not used
0885                       myDetId <<= 5;
0886                       myDetId |= 0;
0887                       // side
0888                       myDetId <<= 2;
0889                       myDetId |= iSide;
0890                       // wheel number
0891                       myDetId <<= 4;
0892                       myDetId |= iWheel;
0893                       // part
0894                       myDetId <<= 2;
0895                       myDetId |= iPart;
0896                       // petal number
0897                       myDetId <<= 4;
0898                       myDetId |= iPetal;
0899                       // ring number
0900                       myDetId <<= 3;
0901                       myDetId |= iRing;
0902                       // module number
0903                       myDetId <<= 3;
0904                       myDetId |= iModule;
0905                       // module type
0906                       myDetId <<= 2;
0907                       myDetId |= iType;
0908                       //
0909                       std::bitset<32> binary_myDetId(myDetId);
0910                       Output << std::endl << std::endl;
0911                       Output << " ******** myDet Id = " << myDetId << " (" << binary_myDetId << ")" << std::endl;
0912                       //
0913                       unsigned int rawid = mapDetIdToGeometricDet[myDetId]->geographicalId().rawId();
0914                       std::bitset<32> binary_detid(rawid);
0915                       GeometricDet::nav_type detNavType = mapDetIdToGeometricDet[myDetId]->navType();
0916                       //
0917                       Output << "            raw Id = " << rawid << " (" << binary_detid << ")"
0918                              << "\t nav type = " << GeometricDet::printNavType(&detNavType.front(), detNavType.size())
0919                              << std::endl;
0920 
0921                       // variables
0922                       fillModuleVariables(mapDetIdToGeometricDet[myDetId], polarRadius, phiRad, z);
0923                       side_z = z;
0924                       wheel_z = z;
0925                       part_z = z;
0926                       switch (iSide) {
0927                           // TEC+
0928                         case 2: {
0929                           //do not change phiRad
0930                           break;
0931                         }
0932                           // TEC-
0933                         case 1: {
0934                           phiRad = changePhiRange_From_ZeroTwoPi_To_MinusPiPlusPi(phiRad);
0935                           break;
0936                         }
0937                         default: {
0938                           // do nothing
0939                         }
0940                       }
0941 
0942                       // petal must be ordered inside [0,2PI) for TEC+, [PI,-PI) for TEC-, take the phi of the central module in a ring with (2n+1) modules
0943                       if ((nModules % 2) && (iModule == (int)(nModules / 2) + (nModules % 2))) {
0944                         switch (iSide) {
0945                             // TEC+
0946                           case 2: {
0947                             petal_phi = phiRad;
0948                             break;
0949                           }
0950                             // TEC-
0951                           case 1: {
0952                             petal_phi = changePhiRange_From_MinusPiPlusPi_To_ZeroTwoPi(phiRad);
0953                             break;
0954                           }
0955                           default: {
0956                             // do nothing
0957                           }
0958                         }
0959                       }
0960 
0961                       ring_R = polarRadius;
0962                       //
0963                       // modules must be ordered inside petals [0,2PI)-->[-PI,PI) if the petal is near phi~0  TEC+ (petal number 1)
0964                       // modules must be ordered inside petals [PI,-PI)-->[2PI,0) if the petal is near phi~PI TEC- (petal number 5)
0965                       switch (iSide) {
0966                           // TEC+
0967                         case 2: {
0968                           if (iPetal == 1) {  // it is the request of the petal at phi = 0, always the first
0969                             module_phi = changePhiRange_From_ZeroTwoPi_To_MinusPiPlusPi(phiRad);
0970                           } else {
0971                             module_phi = phiRad;
0972                           }
0973                           break;
0974                         }
0975                           // TEC-
0976                         case 1: {
0977                           if (iPetal == 5) {  // it is the request of the petal at phi = PI, always the fifth
0978                             module_phi = changePhiRange_From_MinusPiPlusPi_To_MinusTwoPiZero(phiRad);
0979                           } else {
0980                             module_phi = phiRad;
0981                           }
0982                           break;
0983                         }
0984                         default: {
0985                           // do nothing
0986                         }
0987                       }
0988 
0989                       //
0990 
0991                       Output << "\t R = " << polarRadius << "\t phi = " << phiRad << "\t z = " << z << std::endl;
0992 
0993                       // Module Info
0994 
0995                       std::string name = mapDetIdToGeometricDet[myDetId]->name();
0996                       unsigned int theWheel = tTopo->tecWheel(rawid);
0997                       unsigned int theModule = tTopo->tecModule(rawid);
0998                       std::vector<unsigned int> thePetal = tTopo->tecPetalInfo(rawid);
0999                       unsigned int theRing = tTopo->tecRing(rawid);
1000                       std::string side;
1001                       std::string petal;
1002                       side = (tTopo->tecSide(rawid) == 1) ? "-" : "+";
1003                       petal = (thePetal[0] == 1) ? "back" : "front";
1004                       Output << " TEC" << side << "\t"
1005                              << "Wheel " << theWheel << " Petal " << thePetal[1] << " " << petal << " Ring " << theRing
1006                              << "\t"
1007                              << "\t"
1008                              << " module " << theModule << "\t" << name << std::endl;
1009                       //
1010 
1011                       // module: phi check
1012                       Output << "\t # ";
1013                       switch (iSide) {
1014                           // TEC+
1015                         case 2: {
1016                           if ((module_phi - module_phi_previous) < (0 - tolerance_angle)) {
1017                             Output << "\t ERROR phi ordering not respected in module ";
1018                             iERROR++;
1019                           } else {
1020                             Output << "\t OK"
1021                                    << " phi ordering in module ";
1022                             iOK++;
1023                           }
1024                           Output << iModule - 1 << " to " << iModule << " (" << module_phi_previous << " --> "
1025                                  << module_phi << ")" << std::endl;
1026                           break;
1027                         }
1028                           // TEC-
1029                         case 1: {
1030                           if ((module_phi - module_phi_previous) > (0 + tolerance_angle)) {
1031                             Output << "\t ERROR phi ordering not respected in module ";
1032                             iERROR++;
1033                           } else {
1034                             Output << "\t OK"
1035                                    << " phi ordering in module ";
1036                             iOK++;
1037                           }
1038                           Output << iModule - 1 << " to " << iModule << " (" << module_phi_previous << " --> "
1039                                  << module_phi << ")" << std::endl;
1040                           break;
1041                         }
1042                         default: {
1043                           // do nothing
1044                         }
1045                       }
1046                       //
1047                     }  // type loop
1048 
1049                     //
1050                     module_phi_previous = module_phi;
1051                     //
1052                   }  // module loop
1053 
1054                   // ring: R check
1055                   Output << "\t ## ";
1056                   if ((ring_R - ring_R_previous) < (0 + tolerance_space)) {
1057                     Output << "\t ERROR R ordering not respected in petal ";
1058                     iERROR++;
1059                   } else {
1060                     Output << "\t OK"
1061                            << " R ordering in petal ";
1062                     iOK++;
1063                   }
1064                   Output << iPetal << " ring " << iRing - 1 << " to " << iRing << " (" << ring_R_previous << " --> "
1065                          << ring_R << ")" << std::endl;
1066                   //
1067                   switch (iSide) {
1068                       // TEC+
1069                     case 2: {
1070                       ring_R_previous = ring_R;
1071                       module_phi_previous = -M_PI;
1072                       break;
1073                     }
1074                       // TEC-
1075                     case 1: {
1076                       ring_R_previous = ring_R;
1077                       module_phi_previous = 2 * M_PI;
1078                       break;
1079                     }
1080                     default: {
1081                       // do nothing
1082                     }
1083                   }
1084                   //
1085                 }  // ring loop
1086 
1087                 // petal: phi check
1088                 Output << "\t ### ";
1089                 switch (iSide) {
1090                     // TEC+
1091                   case 2: {
1092                     if ((petal_phi - petal_phi_previous) < (0 - tolerance_angle)) {
1093                       Output << "\t ERROR phi ordering not respected in petal ";
1094                       iERROR++;
1095                     } else {
1096                       Output << "\t OK"
1097                              << " phi ordering in petal ";
1098                       iOK++;
1099                     }
1100                     Output << iPetal - 1 << " to " << iPetal << " (" << petal_phi_previous << " --> " << petal_phi
1101                            << ")" << std::endl;
1102                     //
1103                     petal_phi_previous = petal_phi;
1104                     ring_R_previous = 0.0;
1105                     module_phi_previous = -M_PI;
1106                     //
1107                     break;
1108                   }
1109                     // TEC-
1110                   case 1: {
1111                     if ((petal_phi - petal_phi_previous) < (0 - tolerance_angle)) {
1112                       Output << "\t ERROR phi ordering not respected in petal ";
1113                       iERROR++;
1114                     } else {
1115                       Output << "\t OK"
1116                              << " phi ordering in petal ";
1117                       iOK++;
1118                     }
1119                     Output << iPetal - 1 << " to " << iPetal << " (" << petal_phi_previous << " --> " << petal_phi
1120                            << ")" << std::endl;
1121                     //
1122                     petal_phi_previous = petal_phi;
1123                     ring_R_previous = 0.0;
1124                     module_phi_previous = 2 * M_PI;
1125                     //
1126                     break;
1127                   }
1128                   default: {
1129                     // do nothing
1130                   }
1131                 }
1132                 //
1133               }  // petal loop
1134 
1135               // part: |z| check
1136               Output << "\t #### ";
1137               if ((std::abs(part_z) - std::abs(part_z_previous)) < (0 + tolerance_space)) {
1138                 Output << "\t ERROR |z| ordering (front/back) not respected in wheel ";
1139                 iERROR++;
1140               } else {
1141                 Output << "\t OK"
1142                        << " |z| (front/back) ordering in wheel ";
1143                 iOK++;
1144               }
1145               Output << iWheel << " part " << iPart + 1 << " to " << iPart << " (" << part_z_previous << " --> "
1146                      << part_z << ")" << std::endl;
1147               //
1148               switch (iSide) {
1149                   // TEC+
1150                 case 2: {
1151                   part_z_previous = part_z;
1152                   petal_phi_previous = 0.0;
1153                   ring_R_previous = 0.0;
1154                   module_phi_previous = -M_PI;
1155                   break;
1156                 }
1157                   // TEC-
1158                 case 1: {
1159                   part_z_previous = part_z;
1160                   petal_phi_previous = 0.0;
1161                   ring_R_previous = 0.0;
1162                   module_phi_previous = 2 * M_PI;
1163                   break;
1164                 }
1165                 default: {
1166                   // do nothing
1167                 }
1168               }
1169               //
1170             }  // part loop
1171 
1172             // wheel: |z| check
1173             Output << "\t ##### ";
1174             if ((std::abs(wheel_z) - std::abs(wheel_z_previous)) < (0 + tolerance_space)) {
1175               Output << "\t ERROR |z| ordering not respected in wheel ";
1176               iERROR++;
1177             } else {
1178               Output << "\t OK"
1179                      << " |z| ordering in wheel ";
1180               iOK++;
1181             }
1182             Output << iWheel - 1 << " to " << iWheel << " (" << wheel_z_previous << " --> " << wheel_z << ")"
1183                    << std::endl;
1184             //
1185             switch (iSide) {
1186                 // TEC+
1187               case 2: {
1188                 wheel_z_previous = wheel_z;
1189                 part_z_previous = 0.0;
1190                 petal_phi_previous = 0.0;
1191                 ring_R_previous = 0.0;
1192                 module_phi_previous = -M_PI;
1193                 break;
1194               }
1195                 // TEC-
1196               case 1: {
1197                 wheel_z_previous = wheel_z;
1198                 part_z_previous = 0.0;
1199                 petal_phi_previous = 0.0;
1200                 ring_R_previous = 0.0;
1201                 module_phi_previous = 2 * M_PI;
1202                 break;
1203               }
1204               default: {
1205                 // do nothing
1206               }
1207             }
1208             //
1209           }  // wheel loop
1210 
1211           // side: z check
1212           Output << "\t ###### ";
1213           if ((side_z - side_z_previous) > (0 + tolerance_space)) {
1214             Output << "\t ERROR z ordering (side -/+) not respected in TEC side ";
1215             iERROR++;
1216           } else {
1217             Output << "\t OK"
1218                    << " z ordering (side -/+) in TEC side ";
1219             iOK++;
1220           }
1221           Output << iSide + 1 << " to " << iSide << " (" << side_z_previous << " --> " << side_z << ")" << std::endl;
1222           //
1223           switch (iSide) {
1224               // TEC+
1225             case 2: {
1226               side_z_previous = side_z;
1227               wheel_z_previous = 0.0;
1228               part_z_previous = 0.0;
1229               petal_phi_previous = 0.0;
1230               ring_R_previous = 0.0;
1231               module_phi_previous = -M_PI;
1232               break;
1233             }
1234               // TEC-
1235             case 1: {
1236               side_z_previous = side_z;
1237               wheel_z_previous = 0.0;
1238               part_z_previous = 0.0;
1239               petal_phi_previous = 0.0;
1240               ring_R_previous = 0.0;
1241               module_phi_previous = 2 * M_PI;
1242               break;
1243             }
1244             default: {
1245               // do nothing
1246             }
1247           }
1248           //
1249         }  // side loop
1250 
1251         break;
1252       }
1253       default:
1254         Output << " WARNING no Silicon Strip subdetector, I got a " << iSubDetector << std::endl;
1255         ;
1256     }
1257 
1258   }  // subdetector loop
1259 
1260   // summary
1261   unsigned int iChecks = iOK + iERROR;
1262   Output << std::endl << std::endl;
1263   Output << "-------------------------------------" << std::endl;
1264   Output << " Module Numbering Check Summary      " << std::endl;
1265   Output << "-------------------------------------" << std::endl;
1266   Output << " Number of checks:   " << std::setw(6) << iChecks << std::endl;
1267   Output << "               OK:   " << std::setw(6) << iOK << " (" << std::fixed << std::setprecision(2)
1268          << ((double)iOK / (double)iChecks) * 100 << "%) " << std::endl;
1269   Output << "           ERRORS:   " << std::setw(6) << iERROR << " (" << std::fixed << std::setprecision(2)
1270          << ((double)iERROR / (double)iChecks) * 100 << "%) " << std::endl;
1271   Output << "-------------------------------------" << std::endl;
1272 }
1273 
1274 //define this as a plug-in
1275 DEFINE_FWK_MODULE(ModuleNumbering);