1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419
|
/**************************************************************************
* *
* Regina - A Normal Surface Theory Calculator *
* Computational Engine *
* *
* Copyright (c) 1999-2008, Ben Burton *
* For further details contact Ben Burton (bab@debian.org). *
* *
* This program is free software; you can redistribute it and/or *
* modify it under the terms of the GNU General Public License as *
* published by the Free Software Foundation; either version 2 of the *
* License, or (at your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, but *
* WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
* General Public License for more details. *
* *
* You should have received a copy of the GNU General Public *
* License along with this program; if not, write to the Free *
* Software Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, *
* MA 02110-1301, USA. *
* *
**************************************************************************/
/* end stub */
#include <algorithm>
#include "manifold/nsfs.h"
#include "triangulation/ntriangulation.h"
#include "subcomplex/nplugtrisolidtorus.h"
namespace regina {
const int NPlugTriSolidTorus::CHAIN_NONE = 0;
const int NPlugTriSolidTorus::CHAIN_MAJOR = 1;
const int NPlugTriSolidTorus::CHAIN_MINOR = 3;
const int NPlugTriSolidTorus::EQUATOR_MAJOR = 1;
const int NPlugTriSolidTorus::EQUATOR_MINOR = 3;
NPlugTriSolidTorus::~NPlugTriSolidTorus() {
if (core)
delete core;
for (int i = 0; i < 3; i++)
if (chain[i])
delete chain[i];
}
NPlugTriSolidTorus* NPlugTriSolidTorus::clone() const {
NPlugTriSolidTorus* ans = new NPlugTriSolidTorus();
ans->core = core->clone();
for (int i = 0; i < 3; i++) {
if (chain[i])
ans->chain[i] = new NLayeredChain(*chain[i]);
ans->chainType[i] = chainType[i];
}
ans->equatorType = equatorType;
return ans;
}
std::ostream& NPlugTriSolidTorus::writeName(std::ostream& out) const {
long params[3];
int nParams = 0;
int i;
for (i = 0; i < 3; i++)
if (chainType[i] != CHAIN_NONE) {
if (chainType[i] == CHAIN_MAJOR)
params[nParams++] = chain[i]->getIndex();
else
params[nParams++] = -chain[i]->getIndex();
}
std::sort(params, params + nParams);
out << (equatorType == EQUATOR_MAJOR ? "P(" : "P'(");
if (nParams == 0)
return out << "0)";
for (i = 0; i < nParams; i++) {
if (i > 0)
out << ',';
out << params[i];
}
return out << ')';
}
std::ostream& NPlugTriSolidTorus::writeTeXName(std::ostream& out) const {
long params[3];
int nParams = 0;
int i;
for (i = 0; i < 3; i++)
if (chainType[i] != CHAIN_NONE) {
if (chainType[i] == CHAIN_MAJOR)
params[nParams++] = chain[i]->getIndex();
else
params[nParams++] = -chain[i]->getIndex();
}
std::sort(params, params + nParams);
out << (equatorType == EQUATOR_MAJOR ? "P_{" : "P'_{");
if (nParams == 0)
return out << "0}";
for (i = 0; i < nParams; i++) {
if (i > 0)
out << ',';
out << params[i];
}
return out << '}';
}
void NPlugTriSolidTorus::writeTextLong(std::ostream& out) const {
out << "Plugged triangular solid torus: ";
writeName(out);
}
NManifold* NPlugTriSolidTorus::getManifold() const {
NSFSpace* ans = new NSFSpace();
ans->insertFibre(2, -1);
ans->insertFibre(3, 1);
long rot = (equatorType == EQUATOR_MAJOR ? 5 : 4);
for (int i = 0; i < 3; i++)
if (chainType[i] != CHAIN_NONE) {
if (chainType[i] == equatorType)
rot += chain[i]->getIndex();
else
rot -= chain[i]->getIndex();
}
if (rot != 0)
ans->insertFibre(rot, 1);
else {
delete ans;
return 0;
}
ans->reduce();
return ans;
}
NPlugTriSolidTorus* NPlugTriSolidTorus::isPlugTriSolidTorus(
NComponent* comp) {
// TODO: Each triangular solid torus is tested three times since we
// can't call getTetrahedronIndex() from within a component only.
// Basic property checks.
if ((! comp->isClosed()) || (! comp->isOrientable()))
return 0;
if (comp->getNumberOfVertices() > 1)
return 0;
unsigned long nTet = comp->getNumberOfTetrahedra();
if (nTet < 5)
return 0;
// We have a 1-vertex closed orientable component with at least
// 5 tetrahedra.
// Hunt for a core. Make sure we find each triangular solid torus
// just once.
unsigned long tetIndex;
int coreIndex;
NTriSolidTorus* core;
NTetrahedron* coreTet[3];
NEdge* axis[3];
NPerm coreRoles[3];
NTetrahedron* base[2];
NPerm baseRoles[2];
int i, j;
bool error;
NTetrahedron* plugTet[3][2];
NPerm plugRoles[3][2];
NPerm realPlugRoles[2];
NLayeredChain* chain[3];
int chainType[3];
int equatorType = 0;
chain[0] = chain[1] = chain[2] = 0;
for (tetIndex = 0; tetIndex < nTet - 2; tetIndex++)
for (coreIndex = 0; coreIndex < 24; coreIndex++) {
coreRoles[0] = allPermsS4[coreIndex];
if (coreRoles[0][0] > coreRoles[0][3])
continue;
core = NTriSolidTorus::formsTriSolidTorus(
comp->getTetrahedron(tetIndex), coreRoles[0]);
if (! core)
continue;
for (i = 0; i < 3; i++) {
coreTet[i] = core->getTetrahedron(i);
coreRoles[i] = core->getVertexRoles(i);
axis[i] = coreTet[i]->getEdge(
edgeNumber[coreRoles[i][0]][coreRoles[i][3]]);
}
if (axis[0] == axis[1] || axis[1] == axis[2] ||
axis[2] == axis[0]) {
delete core;
continue;
}
// We have the triangular solid torus and we know the three
// axis edges are distinct.
// Hunt for chains.
for (i = 0; i < 3; i++) {
base[0] = coreTet[(i + 1) % 3]->getAdjacentTetrahedron(
coreRoles[(i + 1) % 3][2]);
base[1] = coreTet[(i + 2) % 3]->getAdjacentTetrahedron(
coreRoles[(i + 2) % 3][1]);
if (base[0] != base[1]) {
// No chain.
chainType[i] = CHAIN_NONE;
continue;
}
// Have we layered over the major axis?
baseRoles[0] = coreTet[(i + 1) % 3]->
getAdjacentTetrahedronGluing(coreRoles[(i + 1) % 3][2]) *
coreRoles[(i + 1) % 3] * NPerm(0, 3, 2, 1);
baseRoles[1] = coreTet[(i + 2) % 3]->
getAdjacentTetrahedronGluing(coreRoles[(i + 2) % 3][1]) *
coreRoles[(i + 2) % 3] * NPerm(2, 1, 0, 3);
if (baseRoles[0] == baseRoles[1]) {
chainType[i] = CHAIN_MAJOR;
chain[i] = new NLayeredChain(base[0], baseRoles[0]);
while (chain[i]->extendAbove())
;
continue;
}
// Have we layered over the minor axis?
baseRoles[0] = coreTet[(i + 1) % 3]->
getAdjacentTetrahedronGluing(coreRoles[(i + 1) % 3][2]) *
coreRoles[(i + 1) % 3] * NPerm(3, 0, 2, 1);
baseRoles[1] = coreTet[(i + 2) % 3]->
getAdjacentTetrahedronGluing(coreRoles[(i + 2) % 3][1]) *
coreRoles[(i + 2) % 3] * NPerm(2, 1, 3, 0);
if (baseRoles[0] == baseRoles[1]) {
chainType[i] = CHAIN_MINOR;
chain[i] = new NLayeredChain(base[0], baseRoles[0]);
while (chain[i]->extendAbove())
;
continue;
}
// It's not a chain but it can't be a plug either.
// We'll notice the error because i will be less than 3.
break;
}
// Check whether we broke out of the previous loop with an error.
// Check also whether one of the chains is another in
// reverse, and that we've found the correct number of
// tetrahedra in total.
error = false;
if (i < 3)
error = true;
else if (chain[0] && chain[1] &&
chain[0]->getBottom() == chain[1]->getTop())
error = true;
else if (chain[1] && chain[2] &&
chain[1]->getBottom() == chain[2]->getTop())
error = true;
else if (chain[2] && chain[0] &&
chain[2]->getBottom() == chain[0]->getTop())
error = true;
else if ((chain[0] ? chain[0]->getIndex() : 0) +
(chain[1] ? chain[1]->getIndex() : 0) +
(chain[2] ? chain[2]->getIndex() : 0) +
5 != nTet)
error = true;
if (error) {
for (j = 0; j < 3; j++)
if (chain[j]) {
delete chain[j];
chain[j] = 0;
}
delete core;
continue;
}
// Still hanging in.
// We know there's only 2 tetrahedra left.
// Now we need to check the plug.
error = false;
for (i = 0; i < 3; i++) {
if (chain[i]) {
plugTet[i][0] = chain[i]->getTop()->getAdjacentTetrahedron(
chain[i]->getTopVertexRoles()[3]);
plugTet[i][1] = chain[i]->getTop()->getAdjacentTetrahedron(
chain[i]->getTopVertexRoles()[0]);
plugRoles[i][0] = chain[i]->getTop()->
getAdjacentTetrahedronGluing(chain[i]->
getTopVertexRoles()[3]) *
chain[i]->getTopVertexRoles() *
(chainType[i] == CHAIN_MAJOR ? NPerm(0, 1, 2, 3) :
NPerm(1, 0, 2, 3));
plugRoles[i][1] = chain[i]->getTop()->
getAdjacentTetrahedronGluing(chain[i]->
getTopVertexRoles()[0]) *
chain[i]->getTopVertexRoles() *
(chainType[i] == CHAIN_MAJOR ? NPerm(2, 3, 1, 0) :
NPerm(3, 2, 1, 0));
} else {
plugTet[i][0] = coreTet[(i + 1) % 3]->
getAdjacentTetrahedron(coreRoles[(i + 1) % 3][2]);
plugTet[i][1] = coreTet[(i + 2) % 3]->
getAdjacentTetrahedron(coreRoles[(i + 2) % 3][1]);
plugRoles[i][0] = coreTet[(i + 1) % 3]->
getAdjacentTetrahedronGluing(coreRoles[(i + 1) % 3][2])
* coreRoles[(i + 1) % 3] * NPerm(0, 3, 1, 2);
plugRoles[i][1] = coreTet[(i + 2) % 3]->
getAdjacentTetrahedronGluing(coreRoles[(i + 2) % 3][1])
* coreRoles[(i + 2) % 3] * NPerm(0, 3, 2, 1);
}
}
// Make sure we meet precisely two tetrahedra, three times
// each. Note that this implies that the plug tetrahedra are
// in fact thus far unseen.
for (i = 0; i < 2; i++)
if (plugTet[0][i] != plugTet[1][i] ||
plugTet[1][i] != plugTet[2][i]) {
error = true;
break;
}
// Make sure also that the gluing permutations for the plug
// are correct.
if (! error) {
if (plugRoles[0][0][0] == plugRoles[1][0][0] &&
plugRoles[1][0][0] == plugRoles[2][0][0]) {
// Type EQUATOR_MINOR.
realPlugRoles[0] = plugRoles[0][0] * NPerm(3, 2, 1, 0);
realPlugRoles[1] = plugRoles[0][1] * NPerm(3, 0, 2, 1);
if (realPlugRoles[0] != plugRoles[1][0] *
NPerm(1, 3, 2, 0))
error = true;
else if (realPlugRoles[0] != plugRoles[2][0] *
NPerm(2, 1, 3, 0))
error = true;
else if (realPlugRoles[1] != plugRoles[1][1] *
NPerm(2, 3, 0, 1))
error = true;
else if (realPlugRoles[1] != plugRoles[2][1] *
NPerm(0, 2, 3, 1))
error = true;
else
equatorType = EQUATOR_MINOR;
} else if (plugRoles[0][0][1] == plugRoles[1][0][1] &&
plugRoles[1][0][1] == plugRoles[2][0][1]) {
// Type EQUATOR_MAJOR.
realPlugRoles[0] = plugRoles[0][0] * NPerm(3, 2, 0, 1);
realPlugRoles[1] = plugRoles[0][1] * NPerm(3, 1, 2, 0);
if (realPlugRoles[0] != plugRoles[1][0] *
NPerm(0, 3, 2, 1))
error = true;
else if (realPlugRoles[0] != plugRoles[2][0] *
NPerm(2, 0, 3, 1))
error = true;
else if (realPlugRoles[1] != plugRoles[1][1] *
NPerm(2, 3, 1, 0))
error = true;
else if (realPlugRoles[1] != plugRoles[2][1] *
NPerm(1, 2, 3, 0))
error = true;
else
equatorType = EQUATOR_MAJOR;
} else
error = true;
}
// Finally check the internal face of the plug.
if (! error) {
if (plugTet[0][0]->getAdjacentTetrahedron(realPlugRoles[0][3])
!= plugTet[0][1])
error = true;
else if (plugTet[0][0]->getAdjacentTetrahedronGluing(
realPlugRoles[0][3]) * realPlugRoles[0] !=
realPlugRoles[1])
error = true;
}
if (error) {
for (j = 0; j < 3; j++)
if (chain[j]) {
delete chain[j];
chain[j] = 0;
}
delete core;
continue;
}
// Success!
NPlugTriSolidTorus* plug = new NPlugTriSolidTorus();
plug->core = core;
for (i = 0; i < 3; i++) {
plug->chain[i] = chain[i];
plug->chainType[i] = chainType[i];
}
plug->equatorType = equatorType;
return plug;
}
// Nothing was found.
return 0;
}
} // namespace regina
|