Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-10 02:21:20

0001 ///////////////////////////////////////////////////////////////////////////////
0002 // File: FP420NumberingScheme.cc
0003 // Date: 02.2006
0004 // Description: Numbering scheme for FP420
0005 // Modifications: 08.2008  mside and fside added
0006 ///////////////////////////////////////////////////////////////////////////////
0007 #include "SimG4CMS/FP420/interface/FP420NumberingScheme.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 //
0010 #include <CLHEP/Units/SystemOfUnits.h>
0011 #include "globals.hh"
0012 #include "G4Step.hh"
0013 #include <iostream>
0014 
0015 //#define EDM_ML_DEBUG
0016 
0017 FP420NumberingScheme::FP420NumberingScheme() {
0018   //  sn0=3, pn0=6, rn0=7;
0019 }
0020 
0021 unsigned int FP420NumberingScheme::getUnitID(const G4Step* aStep) const {
0022   unsigned intindex = 0;
0023   const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
0024   int level = (nullptr != touch) ? touch->GetHistoryDepth() + 1 : 0;
0025 
0026 #ifdef EDM_ML_DEBUG
0027   edm::LogVerbatim("FP420") << "FP420NumberingScheme number of levels= " << level;
0028 #endif
0029   if (level > 0) {
0030     int det = 1;
0031     int stationgen = 0;
0032     int zside = 0;
0033     int station = 0;
0034     int plane = 0;
0035     for (int ich = 0; ich < level; ich++) {
0036       int i = level - ich - 1;
0037       const G4String& name = touch->GetVolume(i)->GetName();
0038       int copyno = touch->GetReplicaNumber(i);
0039 
0040       // new and old set up configurations are possible:
0041       if (name == "FP420E") {
0042         det = copyno;
0043       } else if (name == "HPS240E") {
0044         det = copyno + 2;
0045       } else if (name == "FP420Ex1" || name == "HPS240Ex1") {
0046         stationgen = 1;
0047       } else if (name == "FP420Ex3" || name == "HPS240Ex3") {
0048         stationgen = 2;  // was =3
0049       } else if (name == "SISTATION" || name == "HPS240SISTATION") {
0050         station = stationgen;
0051       } else if (name == "SIPLANE" || name == "HPS240SIPLANE") {
0052         plane = copyno;
0053         // SIDETL (or R) can be ether X or Y type in next schemes of readout
0054         //        !!! (=...) zside
0055         //
0056         //      1                  2     <---copyno
0057         //   Front(=2) Empty(=4) Back(=6)     <--SIDETR OR SENSOR2
0058         //      1         2              <---copyno
0059         //   Front(=1) Back(=3) Empty(=5)     <--SIDETL OR SENSOR1
0060         //
0061       } else if (name == "SENSOR2" || name == "HPS240SENSOR2") {
0062         zside = 3 * copyno;  //= 3   6 (copyno=1,2)
0063       } else if (name == "SENSOR1" || name == "HPS240SENSOR1") {
0064         zside = copyno;  //= 1   2 (copyno=1,2)
0065       }
0066 #ifdef EDM_ML_DEBUG
0067       edm::LogVerbatim("FP420") << "FP420NumberingScheme  "
0068                                 << "ich=" << ich << " copyno=" << copyno << " name=" << name;
0069 #endif
0070     }
0071     // det = 1 for +FP420 ,  = 2 for -FP420  / (det-1) = 0,1
0072     // det = 3 for +HPS240 , = 4 for -HPS240 / (det-1) = 2,3
0073     // 0 is as default for every below:
0074     // Z index
0075     // station number 1 - 8   (in reality just 2 ones)
0076     // superplane(superlayer) number  1 - 16 (in reality just 5 ones)
0077 
0078     intindex = packFP420Index(det, zside, station, plane);
0079     //
0080 #ifdef EDM_ML_DEBUG
0081     edm::LogVerbatim("FP420") << "FP420NumberingScheme det=" << det << " zside=" << zside << " station=" << station
0082                               << " plane=" << plane;
0083 #endif
0084   }
0085   return intindex;
0086 }
0087 
0088 unsigned FP420NumberingScheme::packFP420Index(int det, int zside, int station, int plane) {
0089   unsigned int idx = ((det - 1) & 3) << 19;  //bit 19-20 (det-1):0- 3 = 4-->2**2 = 4  -> 4-1  ->((det-1)&3)  2 bit: 0-1
0090   idx += (zside & 7) << 7;                   //bits  7- 9   zside:0- 7 = 8-->2**3 = 8  -> 8-1  ->  (zside&7)  3 bits:0-2
0091   idx += (station & 7) << 4;                 //bits  4- 6 station:0- 7 = 8-->2**3 = 8  -> 8-1  ->(station&7)  3 bits:0-2
0092   idx += (plane & 15);                       //bits  0- 3   plane:0-15 =16-->2**4 =16  ->16-1  -> (plane&15)  4 bits:0-3
0093 
0094 #ifdef EDM_ML_DEBUG
0095   edm::LogVerbatim("FP420") << "FP420 packing: det " << det << " zside  " << zside << " station " << station
0096                             << " plane " << plane << " idx " << idx;
0097 #endif
0098 
0099   return idx;
0100 }
0101 
0102 void FP420NumberingScheme::unpackFP420Index(const unsigned int& idx, int& det, int& zside, int& station, int& plane) {
0103   det = (idx >> 19) & 3;
0104   det += 1;
0105   zside = (idx >> 7) & 7;
0106   station = (idx >> 4) & 7;
0107   plane = idx & 15;
0108   //
0109 
0110 #ifdef EDM_ML_DEBUG
0111   edm::LogVerbatim("FP420") << " FP420unpacking: idx=" << idx << " zside  " << zside << " station " << station
0112                             << " plane " << plane;
0113 #endif
0114   //
0115 }