Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-11 23:27:57

0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002 #include "Geometry/HGCalCommonData/interface/HGCalCell.h"
0003 #include "Geometry/HGCalCommonData/interface/HGCalTypes.h"
0004 #include <vector>
0005 #include <iostream>
0006 
0007 //#define EDM_ML_DEBUG
0008 
0009 HGCalCell::HGCalCell(double waferSize, int32_t nFine, int32_t nCoarse) {
0010   ncell_[0] = nFine;
0011   ncell_[1] = nCoarse;
0012   for (int k = 0; k < 2; ++k) {
0013     cellX_[k] = waferSize / (3 * ncell_[k]);
0014     cellY_[k] = sqrt3By2_ * cellX_[k];
0015   }
0016 #ifdef EDM_ML_DEBUG
0017   edm::LogVerbatim("HGCalGeom") << "HGCalCell initialized with waferSize " << waferSize << " number of cells " << nFine
0018                                 << ":" << nCoarse;
0019 #endif
0020 }
0021 
0022 std::pair<double, double> HGCalCell::cellUV2XY1(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
0023   if (type != 0)
0024     type = 1;
0025   double x(0), y(0);
0026   switch (placementIndex) {
0027     case (HGCalCell::cellPlacementIndex6):
0028       x = (1.5 * (v - u) + 0.5) * cellX_[type];
0029       y = (v + u - 2 * ncell_[type] + 1) * cellY_[type];
0030       break;
0031     case (HGCalCell::cellPlacementIndex7):
0032       x = (1.5 * (v - ncell_[type]) + 1.0) * cellX_[type];
0033       y = (2 * u - v - ncell_[type]) * cellY_[type];
0034       break;
0035     case (HGCalCell::cellPlacementIndex8):
0036       x = (1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
0037       y = -(2 * v - u - ncell_[type] + 1) * cellY_[type];
0038       break;
0039     case (HGCalCell::cellPlacementIndex9):
0040       x = -(1.5 * (v - u) + 0.5) * cellX_[type];
0041       y = -(v + u - 2 * ncell_[type] + 1) * cellY_[type];
0042       break;
0043     case (HGCalCell::cellPlacementIndex10):
0044       x = -(1.5 * (v - ncell_[type]) + 1) * cellX_[type];
0045       y = -(2 * u - v - ncell_[type]) * cellY_[type];
0046       break;
0047     case (HGCalCell::cellPlacementIndex11):
0048       x = -(1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
0049       y = (2 * v - u - ncell_[type] + 1) * cellY_[type];
0050       break;
0051     case (HGCalCell::cellPlacementIndex0):
0052       x = (1.5 * (u - v) - 0.5) * cellX_[type];
0053       y = (v + u - 2 * ncell_[type] + 1) * cellY_[type];
0054       break;
0055     case (HGCalCell::cellPlacementIndex1):
0056       x = -(1.5 * (v - ncell_[type]) + 1.0) * cellX_[type];
0057       y = (2 * u - v - ncell_[type]) * cellY_[type];
0058       break;
0059     case (HGCalCell::cellPlacementIndex2):
0060       x = -(1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
0061       y = -(2 * v - u - ncell_[type] + 1) * cellY_[type];
0062       break;
0063     case (HGCalCell::cellPlacementIndex3):
0064       x = -(1.5 * (u - v) - 0.5) * cellX_[type];
0065       y = -(v + u - 2 * ncell_[type] + 1) * cellY_[type];
0066       break;
0067     case (HGCalCell::cellPlacementIndex4):
0068       x = (1.5 * (v - ncell_[type]) + 1) * cellX_[type];
0069       y = -(2 * u - v - ncell_[type]) * cellY_[type];
0070       break;
0071     default:
0072       x = (1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
0073       y = (2 * v - u - ncell_[type] + 1) * cellY_[type];
0074       break;
0075   }
0076   return std::make_pair(x, y);
0077 }
0078 
0079 std::pair<double, double> HGCalCell::cellUV2XY2(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
0080   if (type != 0)
0081     type = 1;
0082   double x(0), y(0);
0083   if (placementIndex < HGCalCell::cellPlacementExtra) {
0084     double x0 = (1.5 * (u - v) - 0.5) * cellX_[type];
0085     double y0 = (u + v - 2 * ncell_[type] + 1) * cellY_[type];
0086     const std::vector<double> fcos = {1.0, 0.5, -0.5, -1.0, -0.5, 0.5};
0087     const std::vector<double> fsin = {0.0, sqrt3By2_, sqrt3By2_, 0.0, -sqrt3By2_, -sqrt3By2_};
0088     x = x0 * fcos[placementIndex] - y0 * fsin[placementIndex];
0089     y = x0 * fsin[placementIndex] + y0 * fcos[placementIndex];
0090   } else {
0091     double x0 = (1.5 * (v - ncell_[type]) + 1.0) * cellX_[type];
0092     double y0 = (2 * u - v - ncell_[type]) * cellY_[type];
0093     const std::vector<double> fcos = {0.5, 1.0, 0.5, -0.5, -1.0, -0.5};
0094     const std::vector<double> fsin = {sqrt3By2_, 0.0, -sqrt3By2_, -sqrt3By2_, 0.0, sqrt3By2_};
0095     x = x0 * fcos[placementIndex - HGCalCell::cellPlacementExtra] -
0096         y0 * fsin[placementIndex - HGCalCell::cellPlacementExtra];
0097     y = x0 * fsin[placementIndex - HGCalCell::cellPlacementExtra] +
0098         y0 * fcos[placementIndex - HGCalCell::cellPlacementExtra];
0099   }
0100   return std::make_pair(x, y);
0101 }
0102 
0103 std::pair<int, int> HGCalCell::cellUV2Cell(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
0104   if (type != 0)
0105     type = 1;
0106   int cell(0), cellx(0), cellt(HGCalCell::fullCell);
0107   if (placementIndex >= HGCalCell::cellPlacementExtra) {
0108     const std::vector<int> itype0 = {0, 7, 8, 9, 10, 11, 6, 3, 4, 5, 4, 5, 3};
0109     const std::vector<int> itype1 = {0, 0, 1, 2, 3, 4, 5, 0, 1, 2, 0, 1, 2};
0110     const std::vector<int> itype2 = {0, 11, 6, 7, 8, 9, 10, 5, 3, 4, 3, 4, 5};
0111     const std::vector<int> itype3 = {0, 4, 5, 0, 1, 2, 3, 2, 0, 1, 2, 0, 1};
0112     const std::vector<int> itype4 = {0, 9, 10, 11, 6, 7, 8, 4, 5, 3, 5, 3, 4};
0113     const std::vector<int> itype5 = {0, 2, 3, 4, 5, 0, 1, 1, 2, 0, 1, 2, 0};
0114     if (u == 0 && v == 0) {
0115       cellx = 1;
0116       cellt = HGCalCell::cornerCell;
0117     } else if (u == 0 && (v - u) == (ncell_[type] - 1)) {
0118       cellx = 2;
0119       cellt = HGCalCell::cornerCell;
0120     } else if ((v - u) == (ncell_[type] - 1) && v == (2 * ncell_[type] - 1)) {
0121       cellx = 3;
0122       cellt = HGCalCell::cornerCell;
0123     } else if (u == (2 * ncell_[type] - 1) && v == (2 * ncell_[type] - 1)) {
0124       cellx = 4;
0125       cellt = HGCalCell::cornerCell;
0126     } else if (u == (2 * ncell_[type] - 1) && (u - v) == ncell_[type]) {
0127       cellx = 5;
0128       cellt = HGCalCell::cornerCell;
0129     } else if ((u - v) == ncell_[type] && v == 0) {
0130       cellx = 6;
0131       cellt = HGCalCell::cornerCell;
0132     } else if (u == 0) {
0133       cellx = 7;
0134       cellt = HGCalCell::truncatedCell;
0135     } else if ((v - u) == (ncell_[type] - 1)) {
0136       cellx = 10;
0137       cellt = HGCalCell::extendedCell;
0138     } else if (v == (2 * ncell_[type] - 1)) {
0139       cellx = 8;
0140       cellt = HGCalCell::truncatedCell;
0141     } else if (u == (2 * ncell_[type] - 1)) {
0142       cellx = 11;
0143       cellt = HGCalCell::extendedCell;
0144     } else if ((u - v) == ncell_[type]) {
0145       cellx = 9;
0146       cellt = HGCalCell::truncatedCell;
0147     } else if (v == 0) {
0148       cellx = 12;
0149       cellt = HGCalCell::extendedCell;
0150     }
0151     switch (placementIndex) {
0152       case (HGCalCell::cellPlacementIndex6):
0153         cell = itype0[cellx];
0154         break;
0155       case (HGCalCell::cellPlacementIndex7):
0156         cell = itype1[cellx];
0157         break;
0158       case (HGCalCell::cellPlacementIndex8):
0159         cell = itype2[cellx];
0160         break;
0161       case (HGCalCell::cellPlacementIndex9):
0162         cell = itype3[cellx];
0163         break;
0164       case (HGCalCell::cellPlacementIndex10):
0165         cell = itype4[cellx];
0166         break;
0167       default:
0168         cell = itype5[cellx];
0169         break;
0170     }
0171   } else {
0172     const std::vector<int> itype0 = {0, 1, 2, 3, 4, 5, 0, 1, 2, 0, 0, 1, 2};
0173     const std::vector<int> itype1 = {0, 8, 9, 10, 11, 6, 7, 4, 5, 3, 4, 5, 3};
0174     const std::vector<int> itype2 = {0, 3, 4, 5, 0, 1, 2, 2, 0, 1, 1, 2, 0};
0175     const std::vector<int> itype3 = {0, 10, 11, 6, 7, 8, 9, 5, 3, 4, 5, 3, 4};
0176     const std::vector<int> itype4 = {0, 5, 0, 1, 2, 3, 4, 0, 1, 2, 2, 0, 1};
0177     const std::vector<int> itype5 = {0, 6, 7, 8, 9, 10, 11, 3, 4, 5, 3, 4, 5};
0178     if (u == 0 && v == 0) {
0179       cellx = 1;
0180       cellt = HGCalCell::cornerCell;
0181     } else if (v == 0 && (u - v) == (ncell_[type])) {
0182       cellx = 2;
0183       cellt = HGCalCell::cornerCell;
0184     } else if ((u - v) == (ncell_[type]) && u == (2 * ncell_[type] - 1)) {
0185       cellx = 3;
0186       cellt = HGCalCell::cornerCell;
0187     } else if (u == (2 * ncell_[type] - 1) && v == (2 * ncell_[type] - 1)) {
0188       cellx = 4;
0189       cellt = HGCalCell::cornerCell;
0190     } else if (v == (2 * ncell_[type] - 1) && (v - u) == (ncell_[type] - 1)) {
0191       cellx = 5;
0192       cellt = HGCalCell::cornerCell;
0193     } else if ((v - u) == (ncell_[type] - 1) && u == 0) {
0194       cellx = 6;
0195       cellt = HGCalCell::cornerCell;
0196     } else if (v == 0) {
0197       cellx = 10;
0198       cellt = HGCalCell::extendedCell;
0199     } else if ((u - v) == ncell_[type]) {
0200       cellx = 7;
0201       cellt = HGCalCell::truncatedCell;
0202     } else if (u == (2 * ncell_[type] - 1)) {
0203       cellx = 11;
0204       cellt = HGCalCell::extendedCell;
0205     } else if (v == (2 * ncell_[type] - 1)) {
0206       cellx = 8;
0207       cellt = HGCalCell::truncatedCell;
0208     } else if ((v - u) == (ncell_[type] - 1)) {
0209       cellx = 12;
0210       cellt = HGCalCell::extendedCell;
0211     } else if (u == 0) {
0212       cellx = 9;
0213       cellt = HGCalCell::truncatedCell;
0214     }
0215     switch (placementIndex) {
0216       case (HGCalCell::cellPlacementIndex0):
0217         cell = itype0[cellx];
0218         break;
0219       case (HGCalCell::cellPlacementIndex1):
0220         cell = itype1[cellx];
0221         break;
0222       case (HGCalCell::cellPlacementIndex2):
0223         cell = itype2[cellx];
0224         break;
0225       case (HGCalCell::cellPlacementIndex3):
0226         cell = itype3[cellx];
0227         break;
0228       case (HGCalCell::cellPlacementIndex4):
0229         cell = itype4[cellx];
0230         break;
0231       default:
0232         cell = itype5[cellx];
0233         break;
0234     }
0235   }
0236   return std::make_pair(cell, cellt);
0237 }
0238 
0239 int HGCalCell::cellPlacementIndex(int32_t iz, int32_t frontBack, int32_t orient) {
0240   int32_t indx = ((iz * frontBack) > 0) ? orient : (orient + HGCalCell::cellPlacementExtra);
0241   return indx;
0242 }
0243 
0244 std::pair<int32_t, int32_t> HGCalCell::cellOrient(int32_t placementIndex) {
0245   int32_t orient = (placementIndex >= HGCalCell::cellPlacementExtra) ? (placementIndex - HGCalCell::cellPlacementExtra)
0246                                                                      : placementIndex;
0247   int32_t frontBackZside = (placementIndex >= HGCalCell::cellPlacementExtra) ? 1 : -1;
0248   return std::make_pair(orient, frontBackZside);
0249 }
0250 
0251 std::pair<int32_t, int32_t> HGCalCell::cellType(int32_t u, int32_t v, int32_t ncell, int32_t placementIndex) {
0252   int cell(0), cellx(0), cellt(HGCalCell::fullCell);
0253   if (placementIndex >= HGCalCell::cellPlacementExtra) {
0254     const std::vector<int> itype0 = {HGCalCell::centralCell,
0255                                      HGCalCell::bottomCorner,
0256                                      HGCalCell::bottomRightCorner,
0257                                      HGCalCell::topRightCorner,
0258                                      HGCalCell::topCorner,
0259                                      HGCalCell::topLeftCorner,
0260                                      HGCalCell::bottomLeftCorner,
0261                                      HGCalCell::bottomRightEdge,
0262                                      HGCalCell::rightEdge,
0263                                      HGCalCell::topRightEdge,
0264                                      HGCalCell::topLeftEdge,
0265                                      HGCalCell::leftEdge,
0266                                      HGCalCell::bottomLeftEdge};
0267     const std::vector<int> itype1 = {HGCalCell::centralCell,
0268                                      HGCalCell::bottomLeftCorner,
0269                                      HGCalCell::bottomCorner,
0270                                      HGCalCell::bottomRightCorner,
0271                                      HGCalCell::topRightCorner,
0272                                      HGCalCell::topCorner,
0273                                      HGCalCell::topLeftCorner,
0274                                      HGCalCell::bottomLeftEdge,
0275                                      HGCalCell::bottomRightEdge,
0276                                      HGCalCell::rightEdge,
0277                                      HGCalCell::topRightEdge,
0278                                      HGCalCell::topLeftEdge,
0279                                      HGCalCell::leftEdge};
0280     const std::vector<int> itype2 = {HGCalCell::centralCell,
0281                                      HGCalCell::topLeftCorner,
0282                                      HGCalCell::bottomLeftCorner,
0283                                      HGCalCell::bottomCorner,
0284                                      HGCalCell::bottomRightCorner,
0285                                      HGCalCell::topRightCorner,
0286                                      HGCalCell::topCorner,
0287                                      HGCalCell::leftEdge,
0288                                      HGCalCell::bottomLeftEdge,
0289                                      HGCalCell::bottomRightEdge,
0290                                      HGCalCell::rightEdge,
0291                                      HGCalCell::topRightEdge,
0292                                      HGCalCell::topLeftEdge};
0293     const std::vector<int> itype3 = {HGCalCell::centralCell,
0294                                      HGCalCell::topCorner,
0295                                      HGCalCell::topLeftCorner,
0296                                      HGCalCell::bottomLeftCorner,
0297                                      HGCalCell::bottomCorner,
0298                                      HGCalCell::bottomRightCorner,
0299                                      HGCalCell::topRightCorner,
0300                                      HGCalCell::topLeftEdge,
0301                                      HGCalCell::leftEdge,
0302                                      HGCalCell::bottomLeftEdge,
0303                                      HGCalCell::bottomRightEdge,
0304                                      HGCalCell::rightEdge,
0305                                      HGCalCell::topRightEdge};
0306     const std::vector<int> itype4 = {HGCalCell::centralCell,
0307                                      HGCalCell::topRightCorner,
0308                                      HGCalCell::topCorner,
0309                                      HGCalCell::topLeftCorner,
0310                                      HGCalCell::bottomLeftCorner,
0311                                      HGCalCell::bottomCorner,
0312                                      HGCalCell::bottomRightCorner,
0313                                      HGCalCell::topRightEdge,
0314                                      HGCalCell::topLeftEdge,
0315                                      HGCalCell::leftEdge,
0316                                      HGCalCell::bottomLeftEdge,
0317                                      HGCalCell::bottomRightEdge,
0318                                      HGCalCell::rightEdge};
0319     const std::vector<int> itype5 = {HGCalCell::centralCell,
0320                                      HGCalCell::bottomRightCorner,
0321                                      HGCalCell::topRightCorner,
0322                                      HGCalCell::topCorner,
0323                                      HGCalCell::topLeftCorner,
0324                                      HGCalCell::bottomLeftCorner,
0325                                      HGCalCell::bottomCorner,
0326                                      HGCalCell::rightEdge,
0327                                      HGCalCell::topRightEdge,
0328                                      HGCalCell::topLeftEdge,
0329                                      HGCalCell::leftEdge,
0330                                      HGCalCell::bottomLeftEdge,
0331                                      HGCalCell::bottomRightEdge};
0332     if (u == 0 && v == 0) {
0333       cellx = 1;
0334       cellt = HGCalCell::cornerCell;
0335     } else if (u == 0 && (v - u) == (ncell - 1)) {
0336       cellx = 2;
0337       cellt = HGCalCell::cornerCell;
0338     } else if ((v - u) == (ncell - 1) && v == (2 * ncell - 1)) {
0339       cellx = 3;
0340       cellt = HGCalCell::cornerCell;
0341     } else if (u == (2 * ncell - 1) && v == (2 * ncell - 1)) {
0342       cellx = 4;
0343       cellt = HGCalCell::cornerCell;
0344     } else if (u == (2 * ncell - 1) && (u - v) == ncell) {
0345       cellx = 5;
0346       cellt = HGCalCell::cornerCell;
0347     } else if ((u - v) == ncell && v == 0) {
0348       cellx = 6;
0349       cellt = HGCalCell::cornerCell;
0350     } else if (u == 0) {
0351       cellx = 7;
0352       cellt = HGCalCell::truncatedCell;
0353       if (v == 1) {
0354         cellx = 1;
0355         cellt = HGCalCell::truncatedMBCell;
0356       } else if (v == ncell - 2) {
0357         cellx = 2;
0358         cellt = HGCalCell::truncatedMBCell;
0359       }
0360     } else if ((v - u) == (ncell - 1)) {
0361       cellx = 8;
0362       cellt = HGCalCell::extendedCell;
0363       if (v == ncell) {
0364         cellx = 2;
0365         cellt = HGCalCell::extendedMBCell;
0366       } else if (v == 2 * ncell - 2) {
0367         cellx = 3;
0368         cellt = HGCalCell::extendedMBCell;
0369       }
0370     } else if (v == (2 * ncell - 1)) {
0371       cellx = 9;
0372       cellt = HGCalCell::truncatedCell;
0373       if (u == ncell + 1) {
0374         cellx = 3;
0375         cellt = HGCalCell::truncatedMBCell;
0376       } else if (u == 2 * ncell - 2) {
0377         cellx = 4;
0378         cellt = HGCalCell::truncatedMBCell;
0379       }
0380     } else if (u == (2 * ncell - 1)) {
0381       cellx = 10;
0382       cellt = HGCalCell::extendedCell;
0383       if (v == 2 * ncell - 2) {
0384         cellx = 4;
0385         cellt = HGCalCell::extendedMBCell;
0386       } else if (v == ncell + 1) {
0387         cellx = 5;
0388         cellt = HGCalCell::extendedMBCell;
0389       }
0390     } else if ((u - v) == ncell) {
0391       cellx = 11;
0392       cellt = HGCalCell::truncatedCell;
0393       if (u == 2 * ncell - 2) {
0394         cellx = 5;
0395         cellt = HGCalCell::truncatedMBCell;
0396       } else if (u == ncell + 1) {
0397         cellx = 6;
0398         cellt = HGCalCell::truncatedMBCell;
0399       }
0400     } else if (v == 0) {
0401       cellx = 12;
0402       cellt = HGCalCell::extendedCell;
0403       if (u == ncell - 1) {
0404         cellx = 6;
0405         cellt = HGCalCell::extendedMBCell;
0406       } else if (u == 1) {
0407         cellx = 1;
0408         cellt = HGCalCell::extendedMBCell;
0409       }
0410     }
0411     switch (placementIndex) {
0412       case (HGCalCell::cellPlacementIndex6):
0413         cell = itype0[cellx];
0414         break;
0415       case (HGCalCell::cellPlacementIndex7):
0416         cell = itype1[cellx];
0417         break;
0418       case (HGCalCell::cellPlacementIndex8):
0419         cell = itype2[cellx];
0420         break;
0421       case (HGCalCell::cellPlacementIndex9):
0422         cell = itype3[cellx];
0423         break;
0424       case (HGCalCell::cellPlacementIndex10):
0425         cell = itype4[cellx];
0426         break;
0427       default:
0428         cell = itype5[cellx];
0429         break;
0430     }
0431   } else {
0432     const std::vector<int> itype0 = {HGCalCell::centralCell,
0433                                      HGCalCell::bottomCorner,
0434                                      HGCalCell::bottomRightCorner,
0435                                      HGCalCell::topRightCorner,
0436                                      HGCalCell::topCorner,
0437                                      HGCalCell::topLeftCorner,
0438                                      HGCalCell::bottomLeftCorner,
0439                                      HGCalCell::bottomRightEdge,
0440                                      HGCalCell::rightEdge,
0441                                      HGCalCell::topRightEdge,
0442                                      HGCalCell::topLeftEdge,
0443                                      HGCalCell::leftEdge,
0444                                      HGCalCell::bottomLeftEdge};
0445     const std::vector<int> itype1 = {HGCalCell::centralCell,
0446                                      HGCalCell::bottomRightCorner,
0447                                      HGCalCell::topRightCorner,
0448                                      HGCalCell::topCorner,
0449                                      HGCalCell::topLeftCorner,
0450                                      HGCalCell::bottomLeftCorner,
0451                                      HGCalCell::bottomCorner,
0452                                      HGCalCell::rightEdge,
0453                                      HGCalCell::topRightEdge,
0454                                      HGCalCell::topLeftEdge,
0455                                      HGCalCell::leftEdge,
0456                                      HGCalCell::bottomLeftEdge,
0457                                      HGCalCell::bottomRightEdge};
0458     const std::vector<int> itype2 = {HGCalCell::centralCell,
0459                                      HGCalCell::topRightCorner,
0460                                      HGCalCell::topCorner,
0461                                      HGCalCell::topLeftCorner,
0462                                      HGCalCell::bottomLeftCorner,
0463                                      HGCalCell::bottomCorner,
0464                                      HGCalCell::bottomRightCorner,
0465                                      HGCalCell::topRightEdge,
0466                                      HGCalCell::topLeftEdge,
0467                                      HGCalCell::leftEdge,
0468                                      HGCalCell::bottomLeftEdge,
0469                                      HGCalCell::bottomRightEdge,
0470                                      HGCalCell::rightEdge};
0471     const std::vector<int> itype3 = {HGCalCell::centralCell,
0472                                      HGCalCell::topCorner,
0473                                      HGCalCell::topLeftCorner,
0474                                      HGCalCell::bottomLeftCorner,
0475                                      HGCalCell::bottomCorner,
0476                                      HGCalCell::bottomRightCorner,
0477                                      HGCalCell::topRightCorner,
0478                                      HGCalCell::topLeftEdge,
0479                                      HGCalCell::leftEdge,
0480                                      HGCalCell::bottomLeftEdge,
0481                                      HGCalCell::bottomRightEdge,
0482                                      HGCalCell::rightEdge,
0483                                      HGCalCell::topRightEdge};
0484     const std::vector<int> itype4 = {HGCalCell::centralCell,
0485                                      HGCalCell::topLeftCorner,
0486                                      HGCalCell::bottomLeftCorner,
0487                                      HGCalCell::bottomCorner,
0488                                      HGCalCell::bottomRightCorner,
0489                                      HGCalCell::topRightCorner,
0490                                      HGCalCell::topCorner,
0491                                      HGCalCell::leftEdge,
0492                                      HGCalCell::bottomLeftEdge,
0493                                      HGCalCell::bottomRightEdge,
0494                                      HGCalCell::rightEdge,
0495                                      HGCalCell::topRightEdge,
0496                                      HGCalCell::topLeftEdge};
0497     const std::vector<int> itype5 = {HGCalCell::centralCell,
0498                                      HGCalCell::bottomLeftCorner,
0499                                      HGCalCell::bottomCorner,
0500                                      HGCalCell::bottomRightCorner,
0501                                      HGCalCell::topRightCorner,
0502                                      HGCalCell::topCorner,
0503                                      HGCalCell::topLeftCorner,
0504                                      HGCalCell::bottomLeftEdge,
0505                                      HGCalCell::bottomRightEdge,
0506                                      HGCalCell::rightEdge,
0507                                      HGCalCell::topRightEdge,
0508                                      HGCalCell::topLeftEdge,
0509                                      HGCalCell::leftEdge};
0510     if (u == 0 && v == 0) {
0511       cellx = 1;
0512       cellt = HGCalCell::cornerCell;
0513     } else if (v == 0 && (u - v) == (ncell)) {
0514       cellx = 2;
0515       cellt = HGCalCell::cornerCell;
0516     } else if ((u - v) == (ncell) && u == (2 * ncell - 1)) {
0517       cellx = 3;
0518       cellt = HGCalCell::cornerCell;
0519     } else if (u == (2 * ncell - 1) && v == (2 * ncell - 1)) {
0520       cellx = 4;
0521       cellt = HGCalCell::cornerCell;
0522     } else if (v == (2 * ncell - 1) && (v - u) == (ncell - 1)) {
0523       cellx = 5;
0524       cellt = HGCalCell::cornerCell;
0525     } else if ((v - u) == (ncell - 1) && u == 0) {
0526       cellx = 6;
0527       cellt = HGCalCell::cornerCell;
0528     } else if (v == 0) {
0529       cellx = 7;
0530       cellt = HGCalCell::extendedCell;
0531       if (u == 1) {
0532         cellx = 1;
0533         cellt = HGCalCell::extendedMBCell;
0534       } else if (u == ncell - 1) {
0535         cellx = 2;
0536         cellt = HGCalCell::extendedMBCell;
0537       }
0538     } else if ((u - v) == ncell) {
0539       cellx = 8;
0540       cellt = HGCalCell::truncatedCell;
0541       if (u == 2 * ncell - 2) {
0542         cellx = 3;
0543         cellt = HGCalCell::truncatedMBCell;
0544       } else if (u == ncell + 1) {
0545         cellx = 2;
0546         cellt = HGCalCell::truncatedMBCell;
0547       }
0548     } else if (u == (2 * ncell - 1)) {
0549       cellx = 9;
0550       cellt = HGCalCell::extendedCell;
0551       if (v == 2 * ncell - 2) {
0552         cellx = 4;
0553         cellt = HGCalCell::extendedMBCell;
0554       } else if (v == ncell + 1) {
0555         cellx = 3;
0556         cellt = HGCalCell::extendedMBCell;
0557       }
0558     } else if (v == (2 * ncell - 1)) {
0559       cellx = 10;
0560       cellt = HGCalCell::truncatedCell;
0561       if (u == ncell + 1) {
0562         cellx = 5;
0563         cellt = HGCalCell::truncatedMBCell;
0564       } else if (u == 2 * ncell - 2) {
0565         cellx = 4;
0566         cellt = HGCalCell::truncatedMBCell;
0567       }
0568     } else if ((v - u) == (ncell - 1)) {
0569       cellx = 11;
0570       cellt = HGCalCell::extendedCell;
0571       if (v == ncell) {
0572         cellx = 6;
0573         cellt = HGCalCell::extendedMBCell;
0574       } else if (v == 2 * ncell - 2) {
0575         cellx = 5;
0576         cellt = HGCalCell::extendedMBCell;
0577       }
0578     } else if (u == 0) {
0579       cellx = 12;
0580       cellt = HGCalCell::truncatedCell;
0581       if (v == 1) {
0582         cellx = 1;
0583         cellt = HGCalCell::truncatedMBCell;
0584       } else if (v == ncell - 2) {
0585         cellx = 6;
0586         cellt = HGCalCell::truncatedMBCell;
0587       }
0588     }
0589     switch (placementIndex) {
0590       case (HGCalCell::cellPlacementIndex0):
0591         cell = itype0[cellx];
0592         break;
0593       case (HGCalCell::cellPlacementIndex1):
0594         cell = itype1[cellx];
0595         break;
0596       case (HGCalCell::cellPlacementIndex2):
0597         cell = itype2[cellx];
0598         break;
0599       case (HGCalCell::cellPlacementIndex3):
0600         cell = itype3[cellx];
0601         break;
0602       case (HGCalCell::cellPlacementIndex4):
0603         cell = itype4[cellx];
0604         break;
0605       default:
0606         cell = itype5[cellx];
0607         break;
0608     }
0609   }
0610   return std::make_pair(cell, cellt);
0611 }
0612 
0613 std::pair<int32_t, int32_t> HGCalCell::cellType(
0614     int32_t u, int32_t v, int32_t ncell, int32_t placementIndex, int32_t partialType) {
0615   std::pair<int, int> cell = HGCalCell::cellType(u, v, ncell, placementIndex);
0616   int cellx = cell.first;
0617   int cellt = cell.second;
0618   if ((partialType >= HGCalTypes::WaferPartLDOffset) &&
0619       (partialType < (HGCalTypes::WaferPartLDOffset + HGCalTypes::WaferPartLDCount))) {
0620     if ((u == 7 && v == 14) || (u == 7 && v == 0)) {
0621       cellt = HGCalCell::LDPartial0714Cell;
0622       if (u == 7 && v == 0) {
0623         cellx = HGCalCell::leftCell;
0624       } else {
0625         cellx = HGCalCell::rightCell;
0626       }
0627     } else if ((u == 8 && v == 15) || (u == 8 && v == 0)) {
0628       cellt = HGCalCell::LDPartial0815Cell;
0629       if (u == 8 && v == 0) {
0630         cellx = HGCalCell::leftCell;
0631       } else {
0632         cellx = HGCalCell::rightCell;
0633       }
0634     } else if (u == 2 && v == 9) {
0635       cellt = HGCalCell::LDPartial0209Cell;
0636     } else if (u == 0 && v == 7) {
0637       cellt = HGCalCell::LDPartial0007Cell;
0638     } else if (u == 14 && v == 15) {
0639       cellt = HGCalCell::LDPartial1415Cell;
0640     } else if (u == 15 && v == 15) {
0641       cellt = HGCalCell::LDPartial1515Cell;
0642     } else if (u == 7) {
0643       cellt = HGCalCell::extendedCell;
0644       cellx = HGCalCell::topCell;
0645     } else if (u == 8) {
0646       cellt = HGCalCell::truncatedCell;
0647       cellx = HGCalCell::bottomCell;
0648     } else if ((partialType == HGCalTypes::WaferLDLeft) || (partialType == HGCalTypes::WaferLDRight) ||
0649                (partialType == HGCalTypes::WaferLDFive) || (partialType == HGCalTypes::WaferLDThree)) {
0650       if (((u == 15 && v == 11) || (u == 7 && v == 7)) &&
0651           ((partialType == HGCalTypes::WaferLDLeft) || (partialType == HGCalTypes::WaferLDRight))) {
0652         cellt = HGCalCell::halfExtCell;
0653         if (partialType == HGCalTypes::WaferLDLeft) {
0654           cellx = HGCalCell::leftCell;
0655         } else {
0656           cellx = HGCalCell::rightCell;
0657         }
0658       } else if ((u == 7 && v == 11) &&
0659                  ((partialType == HGCalTypes::WaferLDFive) || (partialType == HGCalTypes::WaferLDThree))) {
0660         cellt = HGCalCell::halfExtCell;
0661         if (partialType == HGCalTypes::WaferLDFive) {
0662           cellx = HGCalCell::leftCell;
0663         } else {
0664           cellx = HGCalCell::rightCell;
0665         }
0666       } else if (2 * v - u == 7) {
0667         cellt = HGCalCell::halfCell;
0668         if (partialType == HGCalTypes::WaferLDLeft) {
0669           cellx = HGCalCell::leftCell;
0670         } else if (partialType == HGCalTypes::WaferLDRight) {
0671           cellx = HGCalCell::rightCell;
0672         }
0673       } else if (2 * v - u == 15) {
0674         if (partialType == HGCalTypes::WaferLDFive) {
0675           cellx = HGCalCell::leftCell;
0676         } else if (partialType == HGCalTypes::WaferLDThree) {
0677           cellx = HGCalCell::rightCell;
0678         }
0679         cellt = HGCalCell::halfCell;
0680         std::cout << u << ":" << v << " 1" << std::endl;
0681       }
0682     }
0683   } else if ((partialType >= HGCalTypes::WaferPartHDOffset) &&
0684              (partialType < (HGCalTypes::WaferPartHDOffset + HGCalTypes::WaferPartHDCount))) {
0685     if ((u == 9 && v == 20) || (u == 9 && v == 0)) {
0686       cellt = HGCalCell::HDPartial0920Cell;
0687       if (u == 9 && v == 0) {
0688         cellx = HGCalCell::leftCell;
0689       } else {
0690         cellx = HGCalCell::rightCell;
0691       }
0692     } else if ((u == 10 && v == 21) || (u == 10 && v == 0)) {
0693       cellt = HGCalCell::HDPartial1021Cell;
0694       if (u == 10 && v == 0) {
0695         cellx = HGCalCell::leftCell;
0696       } else {
0697         cellx = HGCalCell::rightCell;
0698       }
0699     } else if (u == 9) {
0700       cellt = HGCalCell::truncatedCell;
0701       cellx = HGCalCell::topCell;
0702     } else if (u == 10) {
0703       cellt = HGCalCell::extendedCell;
0704       cellx = HGCalCell::bottomCell;
0705     } else if ((partialType == HGCalTypes::WaferHDLeft) || (partialType == HGCalTypes::WaferHDRight) ||
0706                (partialType == HGCalTypes::WaferHDFive)) {
0707       if ((u == 10 && v == 7) || (u == 10 && v == 14)) {
0708         cellt = HGCalCell::halfExtCell;
0709         if ((partialType == HGCalTypes::WaferHDLeft) || (partialType == HGCalTypes::WaferHDFive)) {
0710           cellx = HGCalCell::leftCell;
0711         } else {
0712           cellx = HGCalCell::rightCell;
0713         }
0714       } else if ((u == 0 && v == 2) || (u == 0 && v == 9)) {
0715         cellt = HGCalCell::halfTrunCell;
0716         if ((partialType == HGCalTypes::WaferHDLeft) || (partialType == HGCalTypes::WaferHDFive)) {
0717           cellx = HGCalCell::leftCell;
0718         } else {
0719           cellx = HGCalCell::rightCell;
0720         }
0721       } else if (2 * v - u == 4) {
0722         cellt = HGCalCell::halfCell;
0723         if (partialType == HGCalTypes::WaferHDLeft) {
0724           cellx = HGCalCell::leftCell;
0725         }
0726       } else if (2 * v - u == 18) {
0727         cellt = HGCalCell::halfCell;
0728         if (partialType == HGCalTypes::WaferHDFive) {
0729           cellx = HGCalCell::leftCell;
0730         } else if (partialType == HGCalTypes::WaferHDRight) {
0731           cellx = HGCalCell::rightCell;
0732         }
0733       }
0734     }
0735   }
0736   return std::make_pair(cellx, cellt);
0737 }