File indexing completed on 2022-05-28 01:33:52
0001 #include "Geometry/HGCalCommonData/interface/HGCalCell.h"
0002 #include <vector>
0003
0004 HGCalCell::HGCalCell(double waferSize, int32_t nFine, int32_t nCoarse) {
0005 ncell_[0] = nFine;
0006 ncell_[1] = nCoarse;
0007 for (int k = 0; k < 2; ++k) {
0008 cellX_[k] = waferSize / (3 * ncell_[k]);
0009 cellY_[k] = sqrt3By2_ * cellX_[k];
0010 }
0011 }
0012
0013 std::pair<double, double> HGCalCell::cellUV2XY1(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
0014 if (type != 0)
0015 type = 1;
0016 double x(0), y(0);
0017 switch (placementIndex) {
0018 case (HGCalCell::cellPlacementIndex6):
0019 x = (1.5 * (v - u) + 0.5) * cellX_[type];
0020 y = (v + u - 2 * ncell_[type] + 1) * cellY_[type];
0021 break;
0022 case (HGCalCell::cellPlacementIndex7):
0023 x = (1.5 * (v - ncell_[type]) + 1.0) * cellX_[type];
0024 y = (2 * u - v - ncell_[type]) * cellY_[type];
0025 break;
0026 case (HGCalCell::cellPlacementIndex8):
0027 x = (1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
0028 y = -(2 * v - u - ncell_[type] + 1) * cellY_[type];
0029 break;
0030 case (HGCalCell::cellPlacementIndex9):
0031 x = -(1.5 * (v - u) + 0.5) * cellX_[type];
0032 y = -(v + u - 2 * ncell_[type] + 1) * cellY_[type];
0033 break;
0034 case (HGCalCell::cellPlacementIndex10):
0035 x = -(1.5 * (v - ncell_[type]) + 1) * cellX_[type];
0036 y = -(2 * u - v - ncell_[type]) * cellY_[type];
0037 break;
0038 case (HGCalCell::cellPlacementIndex11):
0039 x = -(1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
0040 y = (2 * v - u - ncell_[type] + 1) * cellY_[type];
0041 break;
0042 case (HGCalCell::cellPlacementIndex0):
0043 x = (1.5 * (u - v) - 0.5) * cellX_[type];
0044 y = (v + u - 2 * ncell_[type] + 1) * cellY_[type];
0045 break;
0046 case (HGCalCell::cellPlacementIndex1):
0047 x = -(1.5 * (v - ncell_[type]) + 1.0) * cellX_[type];
0048 y = (2 * u - v - ncell_[type]) * cellY_[type];
0049 break;
0050 case (HGCalCell::cellPlacementIndex2):
0051 x = -(1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
0052 y = -(2 * v - u - ncell_[type] + 1) * cellY_[type];
0053 break;
0054 case (HGCalCell::cellPlacementIndex3):
0055 x = -(1.5 * (u - v) - 0.5) * cellX_[type];
0056 y = -(v + u - 2 * ncell_[type] + 1) * cellY_[type];
0057 break;
0058 case (HGCalCell::cellPlacementIndex4):
0059 x = (1.5 * (v - ncell_[type]) + 1) * cellX_[type];
0060 y = -(2 * u - v - ncell_[type]) * cellY_[type];
0061 break;
0062 default:
0063 x = (1.5 * (u - ncell_[type]) + 0.5) * cellX_[type];
0064 y = (2 * v - u - ncell_[type] + 1) * cellY_[type];
0065 break;
0066 }
0067 return std::make_pair(x, y);
0068 }
0069
0070 std::pair<double, double> HGCalCell::cellUV2XY2(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
0071 if (type != 0)
0072 type = 1;
0073 double x(0), y(0);
0074 if (placementIndex < HGCalCell::cellPlacementExtra) {
0075 double x0 = (1.5 * (u - v) - 0.5) * cellX_[type];
0076 double y0 = (u + v - 2 * ncell_[type] + 1) * cellY_[type];
0077 const std::vector<double> fcos = {1.0, 0.5, -0.5, -1.0, -0.5, 0.5};
0078 const std::vector<double> fsin = {0.0, sqrt3By2_, sqrt3By2_, 0.0, -sqrt3By2_, -sqrt3By2_};
0079 x = x0 * fcos[placementIndex] - y0 * fsin[placementIndex];
0080 y = x0 * fsin[placementIndex] + y0 * fcos[placementIndex];
0081 } else {
0082 double x0 = (1.5 * (v - ncell_[type]) + 1.0) * cellX_[type];
0083 double y0 = (2 * u - v - ncell_[type]) * cellY_[type];
0084 const std::vector<double> fcos = {0.5, 1.0, 0.5, -0.5, -1.0, -0.5};
0085 const std::vector<double> fsin = {sqrt3By2_, 0.0, -sqrt3By2_, -sqrt3By2_, 0.0, sqrt3By2_};
0086 x = x0 * fcos[placementIndex - HGCalCell::cellPlacementExtra] -
0087 y0 * fsin[placementIndex - HGCalCell::cellPlacementExtra];
0088 y = x0 * fsin[placementIndex - HGCalCell::cellPlacementExtra] +
0089 y0 * fcos[placementIndex - HGCalCell::cellPlacementExtra];
0090 }
0091 return std::make_pair(x, y);
0092 }
0093
0094 std::pair<int, int> HGCalCell::cellUV2Cell(int32_t u, int32_t v, int32_t placementIndex, int32_t type) {
0095 if (type != 0)
0096 type = 1;
0097 int cell(0), cellx(0), cellt(HGCalCell::fullCell);
0098 if (placementIndex >= HGCalCell::cellPlacementExtra) {
0099 const std::vector<int> itype0 = {0, 7, 8, 9, 10, 11, 6, 3, 4, 5, 4, 5, 3};
0100 const std::vector<int> itype1 = {0, 0, 1, 2, 3, 4, 5, 0, 1, 2, 0, 1, 2};
0101 const std::vector<int> itype2 = {0, 11, 6, 7, 8, 9, 10, 5, 3, 4, 3, 4, 5};
0102 const std::vector<int> itype3 = {0, 4, 5, 0, 1, 2, 3, 2, 0, 1, 2, 0, 1};
0103 const std::vector<int> itype4 = {0, 9, 10, 11, 6, 7, 8, 4, 5, 3, 5, 3, 4};
0104 const std::vector<int> itype5 = {0, 2, 3, 4, 5, 0, 1, 1, 2, 0, 1, 2, 0};
0105 if (u == 0 && v == 0) {
0106 cellx = 1;
0107 cellt = HGCalCell::cornerCell;
0108 } else if (u == 0 && (v - u) == (ncell_[type] - 1)) {
0109 cellx = 2;
0110 cellt = HGCalCell::cornerCell;
0111 } else if ((v - u) == (ncell_[type] - 1) && v == (2 * ncell_[type] - 1)) {
0112 cellx = 3;
0113 cellt = HGCalCell::cornerCell;
0114 } else if (u == (2 * ncell_[type] - 1) && v == (2 * ncell_[type] - 1)) {
0115 cellx = 4;
0116 cellt = HGCalCell::cornerCell;
0117 } else if (u == (2 * ncell_[type] - 1) && (u - v) == ncell_[type]) {
0118 cellx = 5;
0119 cellt = HGCalCell::cornerCell;
0120 } else if ((u - v) == ncell_[type] && v == 0) {
0121 cellx = 6;
0122 cellt = HGCalCell::cornerCell;
0123 } else if (u == 0) {
0124 cellx = 7;
0125 cellt = HGCalCell::truncatedCell;
0126 } else if ((v - u) == (ncell_[type] - 1)) {
0127 cellx = 10;
0128 cellt = HGCalCell::extendedCell;
0129 } else if (v == (2 * ncell_[type] - 1)) {
0130 cellx = 8;
0131 cellt = HGCalCell::truncatedCell;
0132 } else if (u == (2 * ncell_[type] - 1)) {
0133 cellx = 11;
0134 cellt = HGCalCell::extendedCell;
0135 } else if ((u - v) == ncell_[type]) {
0136 cellx = 9;
0137 cellt = HGCalCell::truncatedCell;
0138 } else if (v == 0) {
0139 cellx = 12;
0140 cellt = HGCalCell::extendedCell;
0141 }
0142 switch (placementIndex) {
0143 case (HGCalCell::cellPlacementIndex6):
0144 cell = itype0[cellx];
0145 break;
0146 case (HGCalCell::cellPlacementIndex7):
0147 cell = itype1[cellx];
0148 break;
0149 case (HGCalCell::cellPlacementIndex8):
0150 cell = itype2[cellx];
0151 break;
0152 case (HGCalCell::cellPlacementIndex9):
0153 cell = itype3[cellx];
0154 break;
0155 case (HGCalCell::cellPlacementIndex10):
0156 cell = itype4[cellx];
0157 break;
0158 default:
0159 cell = itype5[cellx];
0160 break;
0161 }
0162 } else {
0163 const std::vector<int> itype0 = {0, 1, 2, 3, 4, 5, 0, 1, 2, 0, 0, 1, 2};
0164 const std::vector<int> itype1 = {0, 8, 9, 10, 11, 6, 7, 4, 5, 3, 4, 5, 3};
0165 const std::vector<int> itype2 = {0, 3, 4, 5, 0, 1, 2, 2, 0, 1, 1, 2, 0};
0166 const std::vector<int> itype3 = {0, 10, 11, 6, 7, 8, 9, 5, 3, 4, 5, 3, 4};
0167 const std::vector<int> itype4 = {0, 5, 0, 1, 2, 3, 4, 0, 1, 2, 2, 0, 1};
0168 const std::vector<int> itype5 = {0, 6, 7, 8, 9, 10, 11, 3, 4, 5, 3, 4, 5};
0169 if (u == 0 && v == 0) {
0170 cellx = 1;
0171 cellt = HGCalCell::cornerCell;
0172 } else if (v == 0 && (u - v) == (ncell_[type])) {
0173 cellx = 2;
0174 cellt = HGCalCell::cornerCell;
0175 } else if ((u - v) == (ncell_[type]) && u == (2 * ncell_[type] - 1)) {
0176 cellx = 3;
0177 cellt = HGCalCell::cornerCell;
0178 } else if (u == (2 * ncell_[type] - 1) && v == (2 * ncell_[type] - 1)) {
0179 cellx = 4;
0180 cellt = HGCalCell::cornerCell;
0181 } else if (v == (2 * ncell_[type] - 1) && (v - u) == (ncell_[type] - 1)) {
0182 cellx = 5;
0183 cellt = HGCalCell::cornerCell;
0184 } else if ((v - u) == (ncell_[type] - 1) && u == 0) {
0185 cellx = 6;
0186 cellt = HGCalCell::cornerCell;
0187 } else if (v == 0) {
0188 cellx = 10;
0189 cellt = HGCalCell::extendedCell;
0190 } else if ((u - v) == ncell_[type]) {
0191 cellx = 7;
0192 cellt = HGCalCell::truncatedCell;
0193 } else if (u == (2 * ncell_[type] - 1)) {
0194 cellx = 11;
0195 cellt = HGCalCell::extendedCell;
0196 } else if (v == (2 * ncell_[type] - 1)) {
0197 cellx = 8;
0198 cellt = HGCalCell::truncatedCell;
0199 } else if ((v - u) == (ncell_[type] - 1)) {
0200 cellx = 12;
0201 cellt = HGCalCell::extendedCell;
0202 } else if (u == 0) {
0203 cellx = 9;
0204 cellt = HGCalCell::truncatedCell;
0205 }
0206 switch (placementIndex) {
0207 case (HGCalCell::cellPlacementIndex0):
0208 cell = itype0[cellx];
0209 break;
0210 case (HGCalCell::cellPlacementIndex1):
0211 cell = itype1[cellx];
0212 break;
0213 case (HGCalCell::cellPlacementIndex2):
0214 cell = itype2[cellx];
0215 break;
0216 case (HGCalCell::cellPlacementIndex3):
0217 cell = itype3[cellx];
0218 break;
0219 case (HGCalCell::cellPlacementIndex4):
0220 cell = itype4[cellx];
0221 break;
0222 default:
0223 cell = itype5[cellx];
0224 break;
0225 }
0226 }
0227 return std::make_pair(cell, cellt);
0228 }
0229
0230 int HGCalCell::cellPlacementIndex(int32_t iz, int32_t frontBack, int32_t orient) {
0231 int32_t indx = ((iz * frontBack) > 0) ? orient : (orient + HGCalCell::cellPlacementExtra);
0232 return indx;
0233 }
0234
0235 std::pair<int32_t, int32_t> HGCalCell::cellOrient(int32_t placementIndex) {
0236 int32_t orient = (placementIndex >= HGCalCell::cellPlacementExtra) ? (placementIndex - HGCalCell::cellPlacementExtra)
0237 : placementIndex;
0238 int32_t frontBackZside = (placementIndex >= HGCalCell::cellPlacementExtra) ? 1 : -1;
0239 return std::make_pair(orient, frontBackZside);
0240 }
0241
0242 std::pair<int32_t, int32_t> HGCalCell::cellType(int32_t u, int32_t v, int32_t ncell, int32_t placementIndex) {
0243 int cell(0), cellx(0), cellt(HGCalCell::fullCell);
0244 if (placementIndex >= HGCalCell::cellPlacementExtra) {
0245 const std::vector<int> itype0 = {HGCalCell::centralCell,
0246 HGCalCell::bottomCorner,
0247 HGCalCell::bottomRightCorner,
0248 HGCalCell::topRightCorner,
0249 HGCalCell::topCorner,
0250 HGCalCell::topLeftCorner,
0251 HGCalCell::bottomLeftCorner,
0252 HGCalCell::bottomRightEdge,
0253 HGCalCell::rightEdge,
0254 HGCalCell::topRightEdge,
0255 HGCalCell::topLeftEdge,
0256 HGCalCell::leftEdge,
0257 HGCalCell::bottomLeftEdge};
0258 const std::vector<int> itype1 = {HGCalCell::centralCell,
0259 HGCalCell::bottomLeftCorner,
0260 HGCalCell::bottomCorner,
0261 HGCalCell::bottomRightCorner,
0262 HGCalCell::topRightCorner,
0263 HGCalCell::topCorner,
0264 HGCalCell::topLeftCorner,
0265 HGCalCell::bottomLeftEdge,
0266 HGCalCell::bottomRightEdge,
0267 HGCalCell::rightEdge,
0268 HGCalCell::topRightEdge,
0269 HGCalCell::topLeftEdge,
0270 HGCalCell::leftEdge};
0271 const std::vector<int> itype2 = {HGCalCell::centralCell,
0272 HGCalCell::topLeftCorner,
0273 HGCalCell::bottomLeftCorner,
0274 HGCalCell::bottomCorner,
0275 HGCalCell::bottomRightCorner,
0276 HGCalCell::topRightCorner,
0277 HGCalCell::topCorner,
0278 HGCalCell::leftEdge,
0279 HGCalCell::bottomLeftEdge,
0280 HGCalCell::bottomRightEdge,
0281 HGCalCell::rightEdge,
0282 HGCalCell::topRightEdge,
0283 HGCalCell::topLeftEdge};
0284 const std::vector<int> itype3 = {HGCalCell::centralCell,
0285 HGCalCell::topCorner,
0286 HGCalCell::topLeftCorner,
0287 HGCalCell::bottomLeftCorner,
0288 HGCalCell::bottomCorner,
0289 HGCalCell::bottomRightCorner,
0290 HGCalCell::topRightCorner,
0291 HGCalCell::topLeftEdge,
0292 HGCalCell::leftEdge,
0293 HGCalCell::bottomLeftEdge,
0294 HGCalCell::bottomRightEdge,
0295 HGCalCell::rightEdge,
0296 HGCalCell::topRightEdge};
0297 const std::vector<int> itype4 = {HGCalCell::centralCell,
0298 HGCalCell::topRightCorner,
0299 HGCalCell::topCorner,
0300 HGCalCell::topLeftCorner,
0301 HGCalCell::bottomLeftCorner,
0302 HGCalCell::bottomCorner,
0303 HGCalCell::bottomRightCorner,
0304 HGCalCell::topRightEdge,
0305 HGCalCell::topLeftEdge,
0306 HGCalCell::leftEdge,
0307 HGCalCell::bottomLeftEdge,
0308 HGCalCell::bottomRightEdge,
0309 HGCalCell::rightEdge};
0310 const std::vector<int> itype5 = {HGCalCell::centralCell,
0311 HGCalCell::bottomRightCorner,
0312 HGCalCell::topRightCorner,
0313 HGCalCell::topCorner,
0314 HGCalCell::topLeftCorner,
0315 HGCalCell::bottomLeftCorner,
0316 HGCalCell::bottomCorner,
0317 HGCalCell::rightEdge,
0318 HGCalCell::topRightEdge,
0319 HGCalCell::topLeftEdge,
0320 HGCalCell::leftEdge,
0321 HGCalCell::bottomLeftEdge,
0322 HGCalCell::bottomRightEdge};
0323 if (u == 0 && v == 0) {
0324 cellx = 1;
0325 cellt = HGCalCell::cornerCell;
0326 } else if (u == 0 && (v - u) == (ncell - 1)) {
0327 cellx = 2;
0328 cellt = HGCalCell::cornerCell;
0329 } else if ((v - u) == (ncell - 1) && v == (2 * ncell - 1)) {
0330 cellx = 3;
0331 cellt = HGCalCell::cornerCell;
0332 } else if (u == (2 * ncell - 1) && v == (2 * ncell - 1)) {
0333 cellx = 4;
0334 cellt = HGCalCell::cornerCell;
0335 } else if (u == (2 * ncell - 1) && (u - v) == ncell) {
0336 cellx = 5;
0337 cellt = HGCalCell::cornerCell;
0338 } else if ((u - v) == ncell && v == 0) {
0339 cellx = 6;
0340 cellt = HGCalCell::cornerCell;
0341 } else if (u == 0) {
0342 cellx = 7;
0343 cellt = HGCalCell::truncatedCell;
0344 } else if ((v - u) == (ncell - 1)) {
0345 cellx = 10;
0346 cellt = HGCalCell::extendedCell;
0347 } else if (v == (2 * ncell - 1)) {
0348 cellx = 8;
0349 cellt = HGCalCell::truncatedCell;
0350 } else if (u == (2 * ncell - 1)) {
0351 cellx = 11;
0352 cellt = HGCalCell::extendedCell;
0353 } else if ((u - v) == ncell) {
0354 cellx = 9;
0355 cellt = HGCalCell::truncatedCell;
0356 } else if (v == 0) {
0357 cellx = 12;
0358 cellt = HGCalCell::extendedCell;
0359 }
0360 switch (placementIndex) {
0361 case (HGCalCell::cellPlacementIndex6):
0362 cell = itype0[cellx];
0363 break;
0364 case (HGCalCell::cellPlacementIndex7):
0365 cell = itype1[cellx];
0366 break;
0367 case (HGCalCell::cellPlacementIndex8):
0368 cell = itype2[cellx];
0369 break;
0370 case (HGCalCell::cellPlacementIndex9):
0371 cell = itype3[cellx];
0372 break;
0373 case (HGCalCell::cellPlacementIndex10):
0374 cell = itype4[cellx];
0375 break;
0376 default:
0377 cell = itype5[cellx];
0378 break;
0379 }
0380 } else {
0381 const std::vector<int> itype0 = {HGCalCell::centralCell,
0382 HGCalCell::bottomCorner,
0383 HGCalCell::bottomRightCorner,
0384 HGCalCell::topRightCorner,
0385 HGCalCell::topCorner,
0386 HGCalCell::topLeftCorner,
0387 HGCalCell::bottomLeftCorner,
0388 HGCalCell::bottomRightEdge,
0389 HGCalCell::rightEdge,
0390 HGCalCell::topRightEdge,
0391 HGCalCell::topLeftEdge,
0392 HGCalCell::leftEdge,
0393 HGCalCell::bottomLeftEdge};
0394 const std::vector<int> itype1 = {HGCalCell::centralCell,
0395 HGCalCell::bottomRightCorner,
0396 HGCalCell::topRightCorner,
0397 HGCalCell::topCorner,
0398 HGCalCell::topLeftCorner,
0399 HGCalCell::bottomLeftCorner,
0400 HGCalCell::bottomCorner,
0401 HGCalCell::rightEdge,
0402 HGCalCell::topRightEdge,
0403 HGCalCell::topLeftEdge,
0404 HGCalCell::leftEdge,
0405 HGCalCell::bottomLeftEdge,
0406 HGCalCell::bottomRightEdge};
0407 const std::vector<int> itype2 = {HGCalCell::centralCell,
0408 HGCalCell::topRightCorner,
0409 HGCalCell::topCorner,
0410 HGCalCell::topLeftCorner,
0411 HGCalCell::bottomLeftCorner,
0412 HGCalCell::bottomCorner,
0413 HGCalCell::bottomRightCorner,
0414 HGCalCell::topRightEdge,
0415 HGCalCell::topLeftEdge,
0416 HGCalCell::leftEdge,
0417 HGCalCell::bottomLeftEdge,
0418 HGCalCell::bottomRightEdge,
0419 HGCalCell::rightEdge};
0420 const std::vector<int> itype3 = {HGCalCell::centralCell,
0421 HGCalCell::topCorner,
0422 HGCalCell::topLeftCorner,
0423 HGCalCell::bottomLeftCorner,
0424 HGCalCell::bottomCorner,
0425 HGCalCell::bottomRightCorner,
0426 HGCalCell::topRightCorner,
0427 HGCalCell::topLeftEdge,
0428 HGCalCell::leftEdge,
0429 HGCalCell::bottomLeftEdge,
0430 HGCalCell::bottomRightEdge,
0431 HGCalCell::rightEdge,
0432 HGCalCell::topRightEdge};
0433 const std::vector<int> itype4 = {HGCalCell::centralCell,
0434 HGCalCell::topLeftCorner,
0435 HGCalCell::bottomLeftCorner,
0436 HGCalCell::bottomCorner,
0437 HGCalCell::bottomRightCorner,
0438 HGCalCell::topRightCorner,
0439 HGCalCell::topCorner,
0440 HGCalCell::leftEdge,
0441 HGCalCell::bottomLeftEdge,
0442 HGCalCell::bottomRightEdge,
0443 HGCalCell::rightEdge,
0444 HGCalCell::topRightEdge,
0445 HGCalCell::topLeftEdge};
0446 const std::vector<int> itype5 = {HGCalCell::centralCell,
0447 HGCalCell::bottomLeftCorner,
0448 HGCalCell::bottomCorner,
0449 HGCalCell::bottomRightCorner,
0450 HGCalCell::topRightCorner,
0451 HGCalCell::topCorner,
0452 HGCalCell::topLeftCorner,
0453 HGCalCell::bottomLeftEdge,
0454 HGCalCell::bottomRightEdge,
0455 HGCalCell::rightEdge,
0456 HGCalCell::topRightEdge,
0457 HGCalCell::topLeftEdge,
0458 HGCalCell::leftEdge};
0459 if (u == 0 && v == 0) {
0460 cellx = 1;
0461 cellt = HGCalCell::cornerCell;
0462 } else if (v == 0 && (u - v) == (ncell)) {
0463 cellx = 2;
0464 cellt = HGCalCell::cornerCell;
0465 } else if ((u - v) == (ncell) && u == (2 * ncell - 1)) {
0466 cellx = 3;
0467 cellt = HGCalCell::cornerCell;
0468 } else if (u == (2 * ncell - 1) && v == (2 * ncell - 1)) {
0469 cellx = 4;
0470 cellt = HGCalCell::cornerCell;
0471 } else if (v == (2 * ncell - 1) && (v - u) == (ncell - 1)) {
0472 cellx = 5;
0473 cellt = HGCalCell::cornerCell;
0474 } else if ((v - u) == (ncell - 1) && u == 0) {
0475 cellx = 6;
0476 cellt = HGCalCell::cornerCell;
0477 } else if (v == 0) {
0478 cellx = 10;
0479 cellt = HGCalCell::extendedCell;
0480 } else if ((u - v) == ncell) {
0481 cellx = 7;
0482 cellt = HGCalCell::truncatedCell;
0483 } else if (u == (2 * ncell - 1)) {
0484 cellx = 11;
0485 cellt = HGCalCell::extendedCell;
0486 } else if (v == (2 * ncell - 1)) {
0487 cellx = 8;
0488 cellt = HGCalCell::truncatedCell;
0489 } else if ((v - u) == (ncell - 1)) {
0490 cellx = 12;
0491 cellt = HGCalCell::extendedCell;
0492 } else if (u == 0) {
0493 cellx = 9;
0494 cellt = HGCalCell::truncatedCell;
0495 }
0496 switch (placementIndex) {
0497 case (HGCalCell::cellPlacementIndex0):
0498 cell = itype0[cellx];
0499 break;
0500 case (HGCalCell::cellPlacementIndex1):
0501 cell = itype1[cellx];
0502 break;
0503 case (HGCalCell::cellPlacementIndex2):
0504 cell = itype2[cellx];
0505 break;
0506 case (HGCalCell::cellPlacementIndex3):
0507 cell = itype3[cellx];
0508 break;
0509 case (HGCalCell::cellPlacementIndex4):
0510 cell = itype4[cellx];
0511 break;
0512 default:
0513 cell = itype5[cellx];
0514 break;
0515 }
0516 }
0517 return std::make_pair(cell, cellt);
0518 }