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
 420
 421
 422
 423
 424
 425
 426
 427
 428
 429
 430
 431
 432
 433
 434
 435
 436
 437
 438
 439
 440
 441
 442
 443
 444
 445
 446
 447
 448
 449
 450
 451
 452
 453
 454
 455
 456
 457
 458
 459
 460
 461
 462
 463
 464
 465
 466
 467
 468
 469
 470
 471
 472
 473
 474
 475
 476
 477
 478
 479
 480
 481
 482
 483
 484
 485
 486
 487
 488
 489
 490
 491
 492
 493
 494
 495
 496
 497
 498
 499
 500
 501
 502
 503
 504
 505
 506
 507
 508
 509
 510
 511
 512
 513
 514
 515
 516
 517
 518
 519
 520
 521
 522
 523
 524
 525
 526
 527
 528
 529
 530
 531
 532
 533
 534
 535
 536
 537
 538
 539
 540
 541
 542
 543
 544
 545
 546
 547
 548
 549
 550
 551
 552
 553
 554
 555
 556
 557
 558
 559
 560
 561
 562
 563
 564
 565
 566
 567
 568
 569
 570
 571
 572
 573
 574
 575
 576
 577
 578
 579
 580
 581
 582
 583
 584
 585
 586
 587
 588
 589
 590
 591
 592
 593
 594
 595
 596
 597
 598
 599
 600
 601
 602
 603
 604
 605
 606
 607
 608
 609
 610
 611
 612
 613
 614
 615
 616
 617
 618
 619
 620
 621
 622
 623
 624
 625
 626
 627
 628
 629
 630
 631
 632
 633
 634
 635
 636
 637
 638
 639
 640
 641
 642
 643
 644
 645
 646
 647
 648
 649
 650
 651
 652
 653
 654
 655
 656
 657
 658
 659
 660
 661
 662
 663
 664
 665
 666
 667
 668
 669
 670
 671
 672
 673
 674
 675
 676
 677
 678
 679
 680
 681
 682
 683
 684
 685
 686
 687
 688
 689
 690
 691
 692
 693
 694
 695
 696
 697
 698
 699
 700
 701
 702
 703
 704
 705
 706
 707
 708
 709
 710
 711
 712
 713
 714
 715
 716
 717
 718
 719
 720
 721
 722
 723
 724
 725
 726
 727
 728
 729
 730
 731
 732
 733
 734
 735
 736
 737
 738
 739
 740
 741
 742
 743
 744
 745
 746
 747
 748
 749
 750
 751
 752
 753
 754
 755
 756
 757
 758
 759
 760
 761
 762
 763
 764
 765
 766
 767
 768
 769
 770
 771
 772
 773
 774
 775
 776
 777
 778
 779
 780
 781
 782
 783
 784
 785
 786
 787
 788
 789
 790
 791
 792
 793
 794
 795
 796
 797
 798
 799
 800
 801
 802
 803
 804
 805
 806
 807
 808
 809
 810
 811
 812
 813
 814
 815
 816
 817
 818
 819
 820
 821
 822
 823
 824
 825
 826
 827
 828
 829
 830
 831
 832
 833
 834
 835
 836
 837
 838
 839
 840
 841
 842
 843
 844
 845
 846
 847
 848
 849
 850
 851
 852
 853
 854
 855
 856
 857
 858
 859
 860
 861
 862
 863
 864
 865
 866
 867
 868
 869
 870
 871
 872
 873
 874
 875
 876
 877
 878
 879
 880
 881
 882
 883
 884
 885
 886
 887
 888
 889
 890
 891
 892
 893
 894
 895
 896
 897
 898
 899
 900
 901
 902
 903
 904
 905
 906
 907
 908
 909
 910
 911
 912
 913
 914
 915
 916
 917
 918
 919
 920
 921
 922
 923
 924
 925
 926
 927
 928
 929
 930
 931
 932
 933
 934
 935
 936
 937
 938
 939
 940
 941
 942
 943
 944
 945
 946
 947
 948
 949
 950
 951
 952
 953
 954
 955
 956
 957
 958
 959
 960
 961
 962
 963
 964
 965
 966
 967
 968
 969
 970
 971
 972
 973
 974
 975
 976
 977
 978
 979
 980
 981
 982
 983
 984
 985
 986
 987
 988
 989
 990
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
/* *****************************************************************
    MESQUITE -- The Mesh Quality Improvement Toolkit

    Copyright 2004 Sandia Corporation and Argonne National
    Laboratory.  Under the terms of Contract DE-AC04-94AL85000
    with Sandia Corporation, the U.S. Government retains certain
    rights in this software.

    This library is free software; you can redistribute it and/or
    modify it under the terms of the GNU Lesser General Public
    License as published by the Free Software Foundation; either
    version 2.1 of the License, or (at your option) any later version.

    This library 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
    Lesser General Public License for more details.

    You should have received a copy of the GNU Lesser General Public License
    (lgpl.txt) along with this library; if not, write to the Free Software
    Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA

    [email protected], [email protected], [email protected],
    [email protected], [email protected], [email protected]

  ***************************************************************** */
// -*- Mode : c++; tab-width: 3; c-tab-always-indent: t; indent-tabs-mode: nil; c-basic-offset: 3
// -*-

/*!  \File NonSmoothDescent.cpp \brief

  Implements the NonSmoothDescent class member functions.

  \author Lori Freitag
  \date 2002-07-20 */

#include <cstdlib>
#include <cstdio>
#include "NonSmoothDescent.hpp"
#include "MsqTimer.hpp"

#undef NUMERICAL_GRADIENT

namespace MBMesquite
{

const double EPSILON       = 1e-16;
const int MSQ_MAX_OPT_ITER = 20;

enum Rotate
{
    COUNTERCLOCKWISE = 1,
    CLOCKWISE        = 0
};

NonSmoothDescent::NonSmoothDescent( QualityMetric* qm ) : currentQM( qm )<--- Member variable 'NonSmoothDescent::freeVertexIndex' is not initialized in the constructor.<--- Member variable 'NonSmoothDescent::minStepSize' is not initialized in the constructor.<--- Member variable 'NonSmoothDescent::originalValue' is not initialized in the constructor.<--- Member variable 'NonSmoothDescent::optStatus' is not initialized in the constructor.<--- Member variable 'NonSmoothDescent::mSteepest' is not initialized in the constructor.<--- Member variable 'NonSmoothDescent::mAlpha' is not initialized in the constructor.<--- Member variable 'NonSmoothDescent::maxAlpha' is not initialized in the constructor.<--- Member variable 'NonSmoothDescent::pdgInd' is not initialized in the constructor.<--- Member variable 'NonSmoothDescent::mActive' is not initialized in the constructor.<--- Member variable 'NonSmoothDescent::testActive' is not initialized in the constructor.<--- Member variable 'NonSmoothDescent::originalActive' is not initialized in the constructor.
{

    MSQ_DBGOUT( 1 ) << "- Executed NonSmoothDescent::NonSmoothDescent()\n";
}

std::string NonSmoothDescent::get_name() const
{
    return "NonSmoothDescent";
}

PatchSet* NonSmoothDescent::get_patch_set()
{
    return &patchSet;
}

void NonSmoothDescent::initialize( PatchData& /*pd*/, MsqError& /*err*/ )
{
    minStepSize = 1e-6;
    MSQ_DBGOUT( 1 ) << "- Executed NonSmoothDescent::initialize()\n";
}

void NonSmoothDescent::initialize_mesh_iteration( PatchData& /*pd*/, MsqError& /*err*/ ) {}
void NonSmoothDescent::optimize_vertex_positions( PatchData& pd, MsqError& err )
{
    MSQ_FUNCTION_TIMER( "NonSmoothDescent" );

    //  cout << "- Executing NonSmoothDescent::optimize_node_positions()\n";
    /* perform the min max smoothing algorithm */
    MSQ_PRINT( 2 )( "\nInitializing the patch iteration\n" );

    MSQ_PRINT( 3 )( "Number of Vertices: %d\n", (int)pd.num_nodes() );
    MSQ_PRINT( 3 )( "Number of Elements: %d\n", (int)pd.num_elements() );
    // Michael: Note: is this a reliable way to get the dimension?
    // NOTE: Mesquite only does 3-dimensional (including surface) meshes.
    // mDimension = pd.get_mesh()->get_geometric_dimension(err); MSQ_ERRRTN(err);
    // MSQ_PRINT(3)("Spatial Dimension: %d\n",mDimension);
    MSQ_PRINT( 3 )( "Spatial Dimension: 3\n" );

    MSQ_PRINT( 3 )( "Num Free = %d\n", (int)pd.num_free_vertices() );

    MsqFreeVertexIndexIterator free_iter( pd, err );MSQ_ERRRTN( err );
    free_iter.reset();
    free_iter.next();
    freeVertexIndex = free_iter.value();
    MSQ_PRINT( 3 )( "Free Vertex Index = %lu\n", (unsigned long)freeVertexIndex );

    // TODO - need to switch to validity via metric evaluations should
    // be associated with the compute_function somehow
    /* check for an invalid mesh; if it's invalid return and ask the user
       to use untangle */
    if( !this->validity_check( pd, err ) )
    {
        MSQ_PRINT( 1 )( "ERROR: Invalid mesh\n" );
        MSQ_SETERR( err )
        ( "Invalid Mesh: Use untangle to create a valid "
          "triangulation",
          MsqError::INVALID_MESH );
        return;
    }

    /* initialize the optimization data up to numFunctionValues */
    this->init_opt( pd, err );MSQ_ERRRTN( err );
    this->init_max_step_length( pd, err );MSQ_ERRRTN( err );
    MSQ_PRINT( 3 )( "Done initializing optimization\n" );

    /* compute the initial function values */
    // TODO this should return a bool with the validity
    this->compute_function( &pd, originalFunction, err );MSQ_ERRRTN( err );

    // find the initial active set
    this->find_active_set( originalFunction, mActive );

    this->minmax_opt( pd, err );MSQ_ERRRTN( err );
}

void NonSmoothDescent::terminate_mesh_iteration( PatchData& /*pd*/, MsqError& /*err*/ ) {}

void NonSmoothDescent::cleanup()
{
    MSQ_DBGOUT( 1 ) << "- Executing NonSmoothDescent::cleanup()\n";
    MSQ_DBGOUT( 1 ) << "- Done with NonSmoothDescent::cleanup()\n";
}

void NonSmoothDescent::find_plane_points( Direction dir1,
                                          Direction dir2,
                                          const std::vector< Vector3D >& vec,
                                          Vector3D& pt1,
                                          Vector3D& pt2,
                                          Vector3D& /*pt3*/,
                                          Status& status,
                                          MsqError& )
{
    int i;
    int num_min, num_max;
    Rotate rotate   = CLOCKWISE;
    int num_rotated = 0;
    double pt_1, pt_2;
    double min, inv_slope;
    double min_inv_slope = 0.;
    double max;
    double max_inv_slope    = 0;
    double inv_origin_slope = 0;
    const int num_vec       = vec.size();
    const int auto_ind_size = 50;
    int auto_ind[auto_ind_size];
    std::vector< int > heap_ind;
    int* ind;
    if( num_vec <= auto_ind_size )
        ind = auto_ind;
    else
    {
        heap_ind.resize( num_vec );
        ind = &heap_ind[0];
    }

    status = MSQ_CHECK_BOTTOM_UP;
    /* find the minimum points in dir1 starting at -1 */
    num_min = 0;
    ind[0]  = -1;
    ind[1]  = -1;
    ind[2]  = -1;
    min     = 1.0;
    for( i = 0; i < num_vec; i++ )
    {
        if( vec[i][dir1] < min )
        {
            min     = vec[i][dir1];
            ind[0]  = i;
            num_min = 1;
        }
        else if( fabs( vec[i][dir1] - min ) < EPSILON )
        {
            ind[num_min++] = i;
        }
    }
    if( min >= 0 ) status = MSQ_NO_EQUIL;

    if( status != MSQ_NO_EQUIL )
    {
        switch( num_min )
        {
            case 1: /* rotate to find the next point */
                pt1  = vec[ind[0]];
                pt_1 = pt1[dir1];
                pt_2 = pt1[dir2];
                if( pt1[dir2] <= 0 )
                {
                    rotate        = COUNTERCLOCKWISE;
                    max_inv_slope = -HUGE_VAL;
                }
                if( pt1[dir2] > 0 )
                {
                    rotate        = CLOCKWISE;
                    min_inv_slope = HUGE_VAL;
                }
                switch( rotate )
                {
                    case COUNTERCLOCKWISE:
                        for( i = 0; i < num_vec; i++ )
                        {
                            if( i != ind[0] )
                            {
                                inv_slope = ( vec[i][dir2] - pt_2 ) / ( vec[i][dir1] - pt_1 );
                                if( ( inv_slope > max_inv_slope ) && ( fabs( inv_slope - max_inv_slope ) > EPSILON ) )
                                {
                                    ind[1]        = i;
                                    max_inv_slope = inv_slope;
                                    num_rotated   = 1;
                                }
                                else if( fabs( inv_slope - max_inv_slope ) < EPSILON )
                                {
                                    ind[2] = i;
                                    num_rotated++;
                                }
                            }
                        }
                        break;
                    case CLOCKWISE:
                        for( i = 0; i < num_vec; i++ )
                        {
                            if( i != ind[0] )
                            {
                                inv_slope = ( vec[i][dir2] - pt_2 ) / ( vec[i][dir1] - pt_1 );
                                if( ( inv_slope < min_inv_slope ) && ( fabs( inv_slope - max_inv_slope ) > EPSILON ) )
                                {
                                    ind[1]        = i;
                                    min_inv_slope = inv_slope;
                                    num_rotated   = 1;
                                }
                                else if( fabs( inv_slope - min_inv_slope ) < EPSILON )
                                {
                                    ind[2] = i;
                                    num_rotated++;
                                }
                            }
                        }
                }
                switch( num_rotated )
                {
                    case 0:
                        MSQ_PRINT( 3 )( "No points in the rotation ... odd\n" );
                        status = MSQ_HULL_TEST_ERROR;
                        break;
                    case 1:
                        MSQ_PRINT( 3 )( "Found a line in the convex hull\n" );
                        pt2    = vec[ind[1]];
                        status = MSQ_TWO_PT_PLANE;
                        break;
                    default:
                        MSQ_PRINT( 3 )( "Found 2 or more points in the rotation\n" );
                        if( fabs( pt_1 ) > EPSILON ) inv_origin_slope = pt_2 / pt_1;
                        switch( rotate )
                        {
                            case COUNTERCLOCKWISE:
                                if( inv_origin_slope >= max_inv_slope )
                                    status = MSQ_NO_EQUIL;
                                else
                                    status = MSQ_CHECK_TOP_DOWN;
                                break;
                            case CLOCKWISE:
                                if( inv_origin_slope <= min_inv_slope )
                                    status = MSQ_NO_EQUIL;
                                else
                                    status = MSQ_CHECK_TOP_DOWN;
                        }
                }
                break;
            case 2: /* use these two points to define the plane */
                MSQ_PRINT( 3 )( "Found two minimum points to define the plane\n" );
                pt1    = vec[ind[0]];
                pt2    = vec[ind[1]];
                status = MSQ_TWO_PT_PLANE;
                break;
            default: /* check to see if all > 0 */
                MSQ_PRINT( 3 )( "Found 3 or more points in min plane %f\n", min );
                if( vec[ind[0]][dir1] >= 0 )
                    status = MSQ_NO_EQUIL;
                else
                    status = MSQ_CHECK_TOP_DOWN;
        }
    }

    /***************************/
    /*  failed to find any information, checking top/down this coord*/
    /***************************/

    if( status == MSQ_CHECK_TOP_DOWN )
    {
        /* find the maximum points in dir1 starting at 1 */
        num_max = 0;
        ind[0]  = -1;
        ind[1]  = -1;
        ind[2]  = -1;
        max     = -1.0;
        for( i = 0; i < num_vec; i++ )
        {
            if( vec[i][dir1] > max )
            {
                max     = vec[i][dir1];
                ind[0]  = i;
                num_max = 1;
            }
            else if( fabs( vec[i][dir1] - max ) < EPSILON )
            {
                ind[num_max++] = i;
            }
        }
        if( max <= 0 ) status = MSQ_NO_EQUIL;

        if( status != MSQ_NO_EQUIL )
        {
            switch( num_max )
            {
                case 1: /* rotate to find the next point */
                    pt1  = vec[ind[0]];
                    pt_1 = pt1[dir1];
                    pt_2 = pt1[dir2];
                    if( pt1[dir2] < 0 )
                    {
                        rotate        = CLOCKWISE;
                        min_inv_slope = HUGE_VAL;
                    }
                    if( pt1[dir2] >= 0 )
                    {
                        rotate        = COUNTERCLOCKWISE;
                        max_inv_slope = -HUGE_VAL;
                    }
                    switch( rotate )
                    {
                        case COUNTERCLOCKWISE:
                            for( i = 0; i < num_vec; i++ )
                            {
                                if( i != ind[0] )
                                {
                                    inv_slope = ( vec[i][dir2] - pt_2 ) / ( vec[i][dir1] - pt_1 );
                                    if( inv_slope > max_inv_slope )
                                    {
                                        ind[1]        = i;
                                        max_inv_slope = inv_slope;
                                        num_rotated   = 1;
                                    }
                                    else if( fabs( inv_slope - max_inv_slope ) < EPSILON )
                                    {
                                        ind[2] = i;
                                        num_rotated++;
                                    }
                                }
                            }
                            break;
                        case CLOCKWISE:
                            for( i = 0; i < num_vec; i++ )
                            {
                                if( i != ind[0] )
                                {
                                    inv_slope = ( vec[i][dir2] - pt_2 ) / ( vec[i][dir1] - pt_1 );
                                    if( inv_slope < min_inv_slope )
                                    {
                                        ind[1]        = i;
                                        min_inv_slope = inv_slope;
                                        num_rotated   = 1;
                                    }
                                    else if( fabs( inv_slope - min_inv_slope ) < EPSILON )
                                    {
                                        ind[2] = i;
                                        num_rotated++;
                                    }
                                }
                            }
                    }
                    switch( num_rotated )
                    {
                        case 0:
                            MSQ_PRINT( 3 )( "No points in the rotation ... odd\n" );
                            status = MSQ_HULL_TEST_ERROR;
                            break;
                        case 1:
                            MSQ_PRINT( 3 )( "Found a line in the convex hull\n" );
                            pt2    = vec[ind[1]];
                            status = MSQ_TWO_PT_PLANE;
                            break;
                        default:
                            MSQ_PRINT( 3 )( "Found 2 or more points in the rotation\n" );
                            /* check to see if rotation got past origin */
                            inv_origin_slope = pt_2 / pt_1;
                            switch( rotate )
                            {
                                case COUNTERCLOCKWISE:
                                    if( inv_origin_slope >= max_inv_slope )
                                        status = MSQ_NO_EQUIL;
                                    else if( dir1 == 2 )
                                        status = MSQ_CHECK_Y_COORD_DIRECTION;
                                    else if( dir1 == 1 )
                                        status = MSQ_CHECK_X_COORD_DIRECTION;
                                    else
                                        status = MSQ_EQUIL;
                                    break;
                                case CLOCKWISE:
                                    if( inv_origin_slope <= min_inv_slope )
                                        status = MSQ_NO_EQUIL;
                                    else if( dir1 == 2 )
                                        status = MSQ_CHECK_Y_COORD_DIRECTION;
                                    else if( dir1 == 1 )
                                        status = MSQ_CHECK_X_COORD_DIRECTION;
                                    else
                                        status = MSQ_EQUIL;
                            }
                    }
                    break;
                case 2: /* use these two points to define the plane */
                    pt1    = vec[ind[0]];
                    pt2    = vec[ind[1]];
                    status = MSQ_TWO_PT_PLANE;
                    break;
                default: /* check to see if all > 0 */
                    MSQ_PRINT( 3 )( "Found 3 in max plane %f\n", max );
                    if( vec[ind[0]][dir1] <= 0 )
                        status = MSQ_NO_EQUIL;
                    else if( dir1 == 2 )
                        status = MSQ_CHECK_Y_COORD_DIRECTION;
                    else if( dir1 == 1 )
                        status = MSQ_CHECK_X_COORD_DIRECTION;
                    else
                        status = MSQ_EQUIL;
            }
        }
    }
}

void NonSmoothDescent::search_direction( PatchData& /*pd*/, Vector3D& mSearch, MsqError& err )
{
    bool viable;
    double a, b, c, denom;
    std::vector< Vector3D > dir;
    double R0, R1;
    SymmetricMatrix P;
    double x[2];
    double search_mag;

    const int num_active = mActive.active_ind.size();
    ;

    // TODO This might be o.k. actually - i don't see any dependence
    // on the element geometry here... try it and see if it works.
    // if not, try taking all of the gradients in the active set
    // and let the search direction be the average of those.
    //   MSQ_FUNCTION_TIMER( "Search Direction" );

    MSQ_PRINT( 2 )( "\nIn Search Direction\n" );
    this->print_active_set( mActive, mFunction );

    if( num_active == 0 )
    {
        MSQ_SETERR( err )( "No active values in search", MsqError::INVALID_STATE );
        return;
    }

    switch( num_active )
    {
        case 1:
            mSearch   = mGradient[mActive.active_ind[0]];
            mSteepest = mActive.active_ind[0];
            break;
        case 2:
            /* if there are two active points, move in the direction of the
           intersection of the planes.  This is the steepest descent
               direction found by analytically solving the QP */

            /* set up the active gradient directions */
            this->get_active_directions( mGradient, dir, err );MSQ_ERRRTN( err );

            /* form the grammian */
            this->form_grammian( dir, err );MSQ_ERRRTN( err );
            this->form_PD_grammian( err );MSQ_ERRRTN( err );

            denom  = ( mG( 0, 0 ) + mG( 1, 1 ) - 2 * mG( 0, 1 ) );
            viable = true;
            if( fabs( denom ) > EPSILON )
            {
                /* gradients are LI, move along their intersection */
                b = ( mG( 0, 0 ) - mG( 0, 1 ) ) / denom;
                a = 1 - b;
                if( ( b < 0 ) || ( b > 1 ) ) viable = false; /* 0 < b < 1 */
                if( viable )
                {
                    mSearch = a * dir[0] + b * dir[1];
                }
                else
                {
                    /* the gradients are dependent, move along one face */
                    mSearch = dir[0];
                }
            }
            else
            {
                /* the gradients are dependent, move along one face */
                mSearch = dir[0];
            }
            mSteepest = mActive.active_ind[0];

            break;
        default:
            /* as in case 2: solve the QP problem to find the steepest
               descent direction.  This can be done analytically - as
               is done in Gill, Murray and Wright
                 for 3 active points in 3 directions - test PD of G
                 otherwise we know it's SP SD so search edges and faces */

            /* get the active gradient directions */
            this->get_active_directions( mGradient, dir, err );MSQ_ERRRTN( err );

            /* form the entries of the grammian matrix */
            this->form_grammian( dir, err );MSQ_ERRRTN( err );
            this->form_PD_grammian( err );MSQ_ERRRTN( err );

            if( num_active == 3 )
            {
                if( mG.condition3x3() < 1e14 )
                {  // if not singular
                    /* form the entries of P=Z^T G Z where Z = [-1...-1; I ] */
                    this->form_reduced_matrix( P );
                    /* form  the RHS and solve the system for the coeffs */
                    R0      = mG( 0, 0 ) - mG( 1, 0 );
                    R1      = mG( 0, 0 ) - mG( 2, 0 );
                    bool ok = this->solve2x2( P( 0, 0 ), P( 0, 1 ), P( 1, 0 ), P( 1, 1 ), R0, R1, x, err );MSQ_ERRRTN( err );
                    if( ok )
                    {
                        a         = 1 - x[0] - x[1];
                        b         = x[0];
                        c         = x[1];
                        mSearch   = a * dir[0] + b * dir[1] + c * dir[2];
                        mSteepest = mActive.active_ind[0];
                    }
                    else
                    {
                        this->search_edges_faces( &dir[0], mSearch, err );MSQ_ERRRTN( err );
                    }
                }
                else
                {
                    this->search_edges_faces( &dir[0], mSearch, err );MSQ_ERRRTN( err );
                }
            }
            else
            {
                this->search_edges_faces( &dir[0], mSearch, err );MSQ_ERRRTN( err );
            }
    }

    /* if the search direction is essentially zero, equilibrium pt */
    search_mag = mSearch % mSearch;
    MSQ_PRINT( 3 )( "  Search Magnitude %g \n", search_mag );

    if( fabs( search_mag ) < 1E-13 )
        optStatus = MSQ_ZERO_SEARCH;
    else
        mSearch /= std::sqrt( search_mag );
    MSQ_PRINT( 3 )
    ( "  Search Direction %g %g  Steepest %lu\n", mSearch[0], mSearch[1], (unsigned long)mSteepest );
}

void NonSmoothDescent::minmax_opt( PatchData& pd, MsqError& err )
{
    bool equilibriumPt = false;
    Vector3D mSearch( 0, 0, 0 );

    //      int valid;
    MSQ_FUNCTION_TIMER( "Minmax Opt" );
    MSQ_PRINT( 2 )( "In minmax_opt\n" );

    mFunction     = originalFunction;
    originalValue = mActive.true_active_value;

    int iterCount    = 0;
    int optIterCount = 0;

    MSQ_PRINT( 3 )( "Done copying original function to function\n" );

    this->find_active_set( mFunction, mActive );
    prevActiveValues.clear();
    prevActiveValues.push_back( mActive.true_active_value );

    /* check for equilibrium point */
    /* compute the gradient */
    this->compute_gradient( &pd, mGradient, err );MSQ_ERRRTN( err );

    if( mActive.active_ind.size() >= 2 )
    {
        MSQ_PRINT( 3 )( "Testing for an equilibrium point \n" );
        equilibriumPt = this->check_equilibrium( optStatus, err );MSQ_ERRRTN( err );

        if( MSQ_DBG( 2 ) && equilibriumPt ) MSQ_PRINT( 2 )( "Optimization Exiting: An equilibrium point \n" );
    }

    /* terminate if we have found an equilibrium point or if the step is
       too small to be worthwhile continuing */
    while( ( optStatus != MSQ_EQUILIBRIUM ) && ( optStatus != MSQ_STEP_TOO_SMALL ) &&
           ( optStatus != MSQ_IMP_TOO_SMALL ) && ( optStatus != MSQ_FLAT_NO_IMP ) && ( optStatus != MSQ_ZERO_SEARCH ) &&
           ( optStatus != MSQ_MAX_ITER_EXCEEDED ) )
    {

        /* increase the iteration count by one */
        /* smooth_param->iter_count += 1; */
        iterCount += 1;
        optIterCount += 1;
        if( iterCount > MSQ_MAX_OPT_ITER ) optStatus = MSQ_MAX_ITER_EXCEEDED;

        MSQ_PRINT( 3 )( "\nITERATION %d \n", iterCount );

        /* compute the gradient */
        this->compute_gradient( &pd, mGradient, err );MSQ_ERRRTN( err );

        MSQ_PRINT( 3 )( "Computing the search direction \n" );
        this->search_direction( pd, mSearch, err );MSQ_ERRRTN( err );

        /* if there are viable directions to search */
        if( ( optStatus != MSQ_ZERO_SEARCH ) && ( optStatus != MSQ_MAX_ITER_EXCEEDED ) )
        {

            MSQ_PRINT( 3 )( "Computing the projections of the gradients \n" );
            this->get_gradient_projections( mSearch, err );MSQ_ERRRTN( err );

            MSQ_PRINT( 3 )( "Computing the initial step size \n" );
            this->compute_alpha( err );MSQ_ERRRTN( err );

            MSQ_PRINT( 3 )( "Testing whether to accept this step \n" );
            this->step_acceptance( pd, iterCount, mSearch, err );MSQ_ERRRTN( err );
            // MSQ_PRINT(3)("The new free vertex position is %f %f %f\n",
            //  mCoords[freeVertexIndex][0],mCoords[freeVertexIndex][1],mCoords[freeVertexIndex][2]);

            if( MSQ_DBG( 3 ) )
            {
                /* Print the active set */
                this->print_active_set( mActive, mFunction );
            }

            /* check for equilibrium point */
            if( mActive.active_ind.size() >= 2 )
            {
                MSQ_PRINT( 3 )( "Testing for an equilibrium point \n" );
                equilibriumPt = this->check_equilibrium( optStatus, err );MSQ_ERRRTN( err );

                if( MSQ_DBG( 2 ) && equilibriumPt ) MSQ_PRINT( 2 )( "Optimization Exiting: An equilibrium point \n" );
            }

            /* record the values */
            prevActiveValues.push_back( mActive.true_active_value );
        }
        else
        {
            /* decrease the iteration count by one */
            /* smooth_param->iter_count -= 1; */
            iterCount -= 1;
            if( MSQ_DBG( 2 ) )
            {
                MSQ_PRINT( 2 )
                ( "Optimization Exiting: No viable directions; equilibrium point \n" );
                /* Print the old active set */
                this->print_active_set( mActive, mFunction );
            }
        }
    }

    MSQ_PRINT( 2 )( "Checking the validity of the mesh\n" );
    if( !this->validity_check( pd, err ) ) MSQ_PRINT( 2 )( "The final mesh is not valid\n" );MSQ_ERRRTN( err );

    MSQ_PRINT( 2 )( "Number of optimization iterations %d\n", iterCount );

    switch( optStatus )
    {
        default:
            MSQ_PRINT( 2 )( "Optimization Termination OptStatus: Invalid early termination\n" );
            break;
        case MSQ_EQUILIBRIUM:
            MSQ_PRINT( 2 )( "Optimization Termination OptStatus: Equilibrium\n" );
            break;
        case MSQ_STEP_TOO_SMALL:
            MSQ_PRINT( 2 )( "Optimization Termination OptStatus: Step Too Small\n" );
            break;
        case MSQ_IMP_TOO_SMALL:
            MSQ_PRINT( 2 )( "Optimization Termination OptStatus: Improvement Too Small\n" );
            break;
        case MSQ_FLAT_NO_IMP:
            MSQ_PRINT( 2 )( "Optimization Termination OptStatus: Flat No Improvement\n" );
            break;
        case MSQ_ZERO_SEARCH:
            MSQ_PRINT( 2 )( "Optimization Termination OptStatus: Zero Search\n" );
            break;
        case MSQ_MAX_ITER_EXCEEDED:
            MSQ_PRINT( 2 )( "Optimization Termination OptStatus: Max Iter Exceeded\n" );
            break;
    }
}

void NonSmoothDescent::step_acceptance( PatchData& pd, int iterCount, const Vector3D& mSearch, MsqError& err )
{
    const double minAcceptableImprovement = 1e-6;

    //  int        ierr;
    //  int        num_values;
    bool step_done;
    bool valid = true, accept_alpha;
    double estimated_improvement;
    double current_improvement  = HUGE_VAL;
    double previous_improvement = HUGE_VAL;
    double current_percent_diff = HUGE_VAL;<--- Variable 'current_percent_diff' is assigned a value that is never used.
    Vector3D original_point;

    //  MSQ_FUNCTION_TIMER( "Step Acceptance" );
    //  num_values = qmHandles.size();

    optStatus = MSQ_NO_STATUS;

    if( mAlpha < minStepSize )
    {
        optStatus = MSQ_IMP_TOO_SMALL;
        step_done = true;<--- step_done is assigned
        MSQ_PRINT( 3 )( "Alpha starts too small, no improvement\n" );
    }

    const MsqVertex* coords = pd.get_vertex_array( err );

    /* save the original function and active set */
    original_point   = coords[freeVertexIndex];
    originalFunction = mFunction;
    originalActive   = mActive;

    step_done = false;<--- step_done is overwritten
    for( int num_steps = 0; !step_done && num_steps < 100; ++num_steps )
    {
        accept_alpha = false;

        while( !accept_alpha && mAlpha > minStepSize )
        {

            /* make the step */
            pd.move_vertex( -mAlpha * mSearch, freeVertexIndex, err );
            // pd.set_coords_array_element(coords[freeVertexIndex],0,err);

            MSQ_PRINT( 2 )( "search direction %f %f \n", mSearch[0], mSearch[1] );
            MSQ_PRINT( 2 )
            ( "new vertex position %f %f \n", coords[freeVertexIndex][0], coords[freeVertexIndex][1] );

            /* assume alpha is acceptable */
            accept_alpha = true;

            /* never take a step that makes a valid mesh invalid or worsens the quality */
            // TODO Validity check revision -- do the compute function up here
            // and then the rest based on validity
            valid = validity_check( pd, err );MSQ_ERRRTN( err );
            if( valid )
            {
                valid = improvement_check( err );MSQ_ERRRTN( err );
            }
            if( !valid )
            {
                accept_alpha = false;
                pd.move_vertex( mAlpha * mSearch, freeVertexIndex, err );
                // pd.set_coords_array_element(coords[freeVertexIndex],0,err);
                mAlpha = mAlpha / 2;
                MSQ_PRINT( 2 )( "Step not accepted, the new alpha %f\n", mAlpha );

                if( mAlpha < minStepSize )
                {
                    optStatus = MSQ_STEP_TOO_SMALL;
                    step_done = true;
                    MSQ_PRINT( 2 )( "Step too small\n" );
                    /* get back the original point, mFunction, and active set */
                    pd.set_vertex_coordinates( original_point, freeVertexIndex, err );
                    mFunction = originalFunction;
                    mActive   = originalActive;
                }
            }
        }

        if( valid && ( mAlpha > minStepSize ) )
        {
            /* compute the new function and active set */
            this->compute_function( &pd, mFunction, err );MSQ_ERRRTN( err );
            this->find_active_set( mFunction, mActive );

            /* estimate the minimum improvement by taking this step */
            this->get_min_estimate( &estimated_improvement, err );MSQ_ERRRTN( err );
            MSQ_PRINT( 2 )
            ( "The estimated improvement for this step: %f\n", estimated_improvement );

            /* calculate the actual increase */
            current_improvement = mActive.true_active_value - prevActiveValues[iterCount - 1];

            MSQ_PRINT( 3 )( "Actual improvement %f\n", current_improvement );

            /* calculate the percent difference from estimated increase */
            current_percent_diff = fabs( current_improvement - estimated_improvement ) / fabs( estimated_improvement );

            /* determine whether to accept a step */
            if( ( fabs( previous_improvement ) > fabs( current_improvement ) ) && ( previous_improvement < 0 ) )
            {
                /* accept the previous step - it was better */
                MSQ_PRINT( 2 )( "Accepting the previous step\n" );

                /* subtract alpha in again (previous step) */
                pd.move_vertex( -mAlpha * mSearch, freeVertexIndex, err );
                // pd.set_coords_array_element(coords[freeVertexIndex],0,err);

                /* does this make an invalid mesh valid? */
                // TODO Validity check revisison
                valid = validity_check( pd, err );MSQ_ERRRTN( err );
                if( valid ) valid = improvement_check( err );MSQ_ERRRTN( err );

                /* copy test function and active set */
                mFunction = testFunction;
                mActive   = testActive;

                optStatus = MSQ_STEP_ACCEPTED;
                step_done = true;

                /* check to see that we're still making good improvements */
                if( fabs( previous_improvement ) < minAcceptableImprovement )
                {
                    optStatus = MSQ_IMP_TOO_SMALL;
                    step_done = true;
                    MSQ_PRINT( 2 )( "Optimization Exiting: Improvement too small\n" );
                }
            }
            else if( ( ( fabs( current_improvement ) > fabs( estimated_improvement ) ) ||
                       ( current_percent_diff < .1 ) ) &&
                     ( current_improvement < 0 ) )
            {
                /* accept this step, exceeded estimated increase or was close */
                optStatus = MSQ_STEP_ACCEPTED;
                step_done = true;

                /* check to see that we're still making good improvements */
                if( fabs( current_improvement ) < minAcceptableImprovement )
                {
                    MSQ_PRINT( 2 )( "Optimization Exiting: Improvement too small\n" );
                    optStatus = MSQ_IMP_TOO_SMALL;
                    step_done = true;
                }
            }
            else if( ( current_improvement > 0 ) && ( previous_improvement > 0 ) &&
                     ( fabs( current_improvement ) < minAcceptableImprovement ) &&
                     ( fabs( previous_improvement ) < minAcceptableImprovement ) )
            {

                /* we are making no progress, quit */
                optStatus = MSQ_FLAT_NO_IMP;
                step_done = true;
                MSQ_PRINT( 2 )( "Opimization Exiting: Flat no improvement\n" );

                /* get back the original point, function, and active set */
                pd.set_vertex_coordinates( original_point, freeVertexIndex, err );MSQ_ERRRTN( err );
                mFunction = originalFunction;
                mActive   = originalActive;
            }
            else
            {
                /* halve alpha and try again */
                /* add out the old step */
                pd.move_vertex( mAlpha * mSearch, freeVertexIndex, err );
                // pd.set_coords_array_element(coords[freeVertexIndex],0,err);

                /* halve step size */
                mAlpha = mAlpha / 2;
                MSQ_PRINT( 3 )( "Step not accepted, the new alpha %f\n", mAlpha );

                if( mAlpha < minStepSize )
                {
                    /* get back the original point, function, and active set */
                    MSQ_PRINT( 2 )( "Optimization Exiting: Step too small\n" );
                    pd.set_vertex_coordinates( original_point, freeVertexIndex, err );MSQ_ERRRTN( err );
                    mFunction = originalFunction;
                    mActive   = originalActive;
                    optStatus = MSQ_STEP_TOO_SMALL;
                    step_done = true;
                }
                else
                {
                    testFunction         = mFunction;
                    testActive           = mActive;
                    previous_improvement = current_improvement;
                }
            }
        }
    }
    if( current_improvement > 0 && optStatus == MSQ_STEP_ACCEPTED )
    {
        MSQ_PRINT( 2 )( "Accepted a negative step %f \n", current_improvement );
    }
}

bool NonSmoothDescent::compute_function( PatchData* patch_data, std::vector< double >& func, MsqError& err )
{
    // NEED 1.0/FUNCTION WHICH IS ONLY TRUE OF CONDITION NUMBER

    func.resize( qmHandles.size() );
    bool valid_bool = true;
    for( size_t i = 0; i < qmHandles.size(); i++ )
    {
        valid_bool = valid_bool && currentQM->evaluate( *patch_data, qmHandles[i], func[i], err );
        MSQ_ERRZERO( err );
    }

    return valid_bool;
}

bool NonSmoothDescent::compute_gradient( PatchData* patch_data, std::vector< Vector3D >& gradient_out, MsqError& err )
{
    MSQ_DBGOUT( 2 ) << "Computing Gradient\n";

    bool valid = true;

#ifdef NUMERICAL_GRADIENT

    const double delta = 10e-6;
    std::vector< double > func( qmHandles.size() ), fdelta( qmHandles.size() );

    valid = this->compute_function( patch_data, func, err );
    MSQ_ERRZERO( err );

    /* gradient in the x, y, z direction */

    for( int j = 0; j < 3; j++ )
    {

        // perturb the coordinates of the free vertex in the j direction by delta
        Vector3D delta_3( 0, 0, 0 );
        Vector3D orig_pos = patch_data->vertex_by_index( freeVertexIndex );
        delta_3[j]        = delta;
        patch_data->move_vertex( delta_3, freeVertexIndex, err );

        // compute the function at the perturbed point location
        valid = valid && this->compute_function( patch_data, fdelta, err );
        MSQ_ERRZERO( err );

        // compute the numerical gradient
        for( size_t i = 0; i < qmHandles.size(); i++ )
        {
            mGradient[i][j] = ( fdelta[i] - func[i] ) / delta;
            // MSQ_DEBUG_ACTION(3,{fprintf(stdout,"  Gradient
            // value[%d][%d]=%g\n",i,j,mGradient[i][j]);});
        }

        // put the coordinates back where they belong
        patch_data->set_vertex_coordinates( orig_pos, freeVertexIndex, err );
    }

#else
    double value;
    gradient_out.resize( qmHandles.size() );
    for( size_t i = 0; i < qmHandles.size(); i++ )
    {
        valid = valid && currentQM->evaluate_with_gradient( *patch_data, qmHandles[i], value, tmpIdx, tmpGrad, err );
        MSQ_ERRZERO( err );
        assert( tmpIdx.size() == 1 && tmpIdx[0] == freeVertexIndex );
        gradient_out[i] = tmpGrad[0];
    }
#endif

    return valid;
}

void NonSmoothDescent::find_active_set( const std::vector< double >& function, ActiveSet& active_set )
{
    // local parameter initialization
    const double activeEpsilon = .3e-4;
    //  activeEpsilon = .3e-8;

    double function_val;
    double active_value0;
    double temp;

    //    FUNCTION_TIMER_START("Find Active Set");
    MSQ_DBGOUT( 2 ) << "\nFinding the active set\n";

    /* the first function value is our initial active value */
    active_set.set( 0 );
    active_set.true_active_value = function[0];
    //    MSQ_DEBUG_ACTION(3,{fprintf(stdout,"  function value[0]: %g\n",function[0]);});

    /* first sort out the active set...
       all vals within active_epsilon of largest val */

    for( size_t i = 1; i < qmHandles.size(); i++ )
    {
        function_val = function[i];
        if( active_set.true_active_value < function_val ) active_set.true_active_value = function_val;
        active_value0 = function[active_set.active_ind[0]];
        temp          = fabs( function_val - active_value0 );
        //        MSQ_DEBUG_ACTION(3,{fprintf(stdout,"  function value[%d]: %g\n",i,function[i]);});
        if( function_val > active_value0 )
        {  // seek max_i function[i]
            if( temp >= activeEpsilon )
            {
                active_set.set( i );  // new max
            }
            else
            {
                active_set.add( i, fabs( function_val - active_value0 ) < EPSILON );
            }
        }
        else
        {
            if( temp < activeEpsilon )
            {
                active_set.add( i, fabs( function_val - active_value0 ) < EPSILON );
            }
        }
    }
}

bool NonSmoothDescent::validity_check( PatchData& pd, MsqError& err )

{
    //  FUNCTION_TIMER_START("Validity Check");

    // ONLY FOR SIMPLICIAL MESHES - THERE SHOULD BE A VALIDITY CHECKER ASSOCIATED
    // WITH MSQ ELEMENTS

    /* check that the simplicial mesh is still valid, based on right handedness.
         Returns a 1 or a 0 */

    // TODO as a first step we can switch this over to the function
    // evaluation and use the rest of the code as is
    bool valid  = true;
    double dEps = 1.e-13;

    double x1, x2, x3, x4, y1, y2, y3, y4, z1, z2, z3, z4;
    const MsqVertex* coords = pd.get_vertex_array( err );

    for( size_t i = 0; i < pd.num_elements(); i++ )
    {
        const size_t* conn = pd.element_by_index( i ).get_vertex_index_array();
        coords[conn[0]].get_coordinates( x1, y1, z1 );
        coords[conn[1]].get_coordinates( x2, y2, z2 );
        coords[conn[2]].get_coordinates( x3, y3, z3 );
        coords[conn[3]].get_coordinates( x4, y4, z4 );

        double dDX2 = x2 - x1;
        double dDX3 = x3 - x1;
        double dDX4 = x4 - x1;

        double dDY2 = y2 - y1;
        double dDY3 = y3 - y1;
        double dDY4 = y4 - y1;

        double dDZ2 = z2 - z1;
        double dDZ3 = z3 - z1;
        double dDZ4 = z4 - z1;

        /* dDet is proportional to the cell volume */
        double dDet = dDX2 * dDY3 * dDZ4 + dDX3 * dDY4 * dDZ2 + dDX4 * dDY2 * dDZ3 - dDZ2 * dDY3 * dDX4 -
                      dDZ3 * dDY4 * dDX2 - dDZ4 * dDY2 * dDX3;

        /* Compute a length scale based on edge lengths. */
        double dScale = ( sqrt( ( x1 - x2 ) * ( x1 - x2 ) + ( y1 - y2 ) * ( y1 - y2 ) + ( z1 - z2 ) * ( z1 - z2 ) ) +
                          sqrt( ( x1 - x3 ) * ( x1 - x3 ) + ( y1 - y3 ) * ( y1 - y3 ) + ( z1 - z3 ) * ( z1 - z3 ) ) +
                          sqrt( ( x1 - x4 ) * ( x1 - x4 ) + ( y1 - y4 ) * ( y1 - y4 ) + ( z1 - z4 ) * ( z1 - z4 ) ) +
                          sqrt( ( x2 - x3 ) * ( x2 - x3 ) + ( y2 - y3 ) * ( y2 - y3 ) + ( z2 - z3 ) * ( z2 - z3 ) ) +
                          sqrt( ( x2 - x4 ) * ( x2 - x4 ) + ( y2 - y4 ) * ( y2 - y4 ) + ( z2 - z4 ) * ( z2 - z4 ) ) +
                          sqrt( ( x3 - x4 ) * ( x3 - x4 ) + ( y3 - y4 ) * ( y3 - y4 ) + ( z3 - z4 ) * ( z3 - z4 ) ) ) /
                        6.;

        /* Use the length scale to get a better idea if the tet is flat or
           just really small. */
        if( fabs( dScale ) < EPSILON )
        {
            return false;
        }
        else
        {
            dDet /= ( dScale * dScale * dScale );
        }

        if( dDet > dEps )
        {
            valid = true;
        }
        else if( dDet < -dEps )
        {
            return false;
        }
        else
        {
            return false;
        }
    }  // end for i=1,numElements

    //  MSQ_DEBUG_ACTION(2,{fprintf(stdout,"Mesh Validity is: %d \n",valid);});

    //  FUNCTION_TIMER_END();
    return ( valid );
}

void NonSmoothDescent::get_active_directions( const std::vector< Vector3D >& pGradient,
                                              std::vector< Vector3D >& dir,
                                              MsqError& /*err*/ )
{
    const size_t num_active = mActive.active_ind.size();
    dir.resize( num_active );
    for( size_t i = 0; i < num_active; i++ )
    {
        dir[i] = pGradient[mActive.active_ind[i]];
    }
}

bool NonSmoothDescent::check_vector_dots( const std::vector< Vector3D >& vec,
                                          const Vector3D& normal,
                                          MsqError& /*err*/ )
{
    double test_dot = vec[0] % normal;
    unsigned ind    = 1;
    while( ( fabs( test_dot ) < EPSILON ) && ( ind < vec.size() ) )
    {
        test_dot = vec[ind] % normal;
        ind++;
    }

    for( unsigned i = ind; i < vec.size(); i++ )
    {
        double dot = vec[i] % normal;
        if( ( ( dot > 0 && test_dot < 0 ) || ( dot < 0 && test_dot > 0 ) ) && ( fabs( dot ) > EPSILON ) )
        {
            return true;
        }
    }
    return false;
}

bool NonSmoothDescent::convex_hull_test( const std::vector< Vector3D >& vec, MsqError& err )
{
    //    int ierr;
    bool equil         = false;
    Direction dir_done = MSQ_ZDIR;
    Status status      = MSQ_CHECK_Z_COORD_DIRECTION;
    Vector3D pt1, pt2, pt3, normal;

    /* tries to determine equilibrium for the 3D case */

    if( vec.size() <= 2 ) status = MSQ_NO_EQUIL;

    while( ( status != MSQ_EQUIL ) && ( status != MSQ_NO_EQUIL ) && ( status != MSQ_HULL_TEST_ERROR ) )
    {
        if( status == MSQ_CHECK_Z_COORD_DIRECTION )
        {
            this->find_plane_points( MSQ_ZDIR, MSQ_YDIR, vec, pt1, pt2, pt3, status, err );
            dir_done = MSQ_ZDIR;
        }
        if( status == MSQ_CHECK_Y_COORD_DIRECTION )
        {
            this->find_plane_points( MSQ_YDIR, MSQ_XDIR, vec, pt1, pt2, pt3, status, err );
            dir_done = MSQ_YDIR;
        }
        if( status == MSQ_CHECK_X_COORD_DIRECTION )
        {
            this->find_plane_points( MSQ_XDIR, MSQ_ZDIR, vec, pt1, pt2, pt3, status, err );
            dir_done = MSQ_XDIR;
        }
        if( status == MSQ_TWO_PT_PLANE )
        {
            pt3 = Vector3D( 0, 0, 0 );
        }
        if( ( status == MSQ_TWO_PT_PLANE ) || ( status == MSQ_THREE_PT_PLANE ) )
        {
            this->find_plane_normal( pt1, pt2, pt3, normal, err );
            equil = this->check_vector_dots( vec, normal, err );
            if( equil )<--- Assuming that condition 'equil' is not redundant
            {
                switch( dir_done )
                {
                    case MSQ_ZDIR:
                        equil  = false;
                        status = MSQ_CHECK_Y_COORD_DIRECTION;
                        break;
                    case MSQ_YDIR:
                        equil  = false;
                        status = MSQ_CHECK_X_COORD_DIRECTION;
                        break;
                    case MSQ_XDIR:
                        equil  = true;
                        status = MSQ_EQUIL;
                }
            }
            else if( equil == 0 )<--- Condition 'equil==0' is always true
            {
                status = MSQ_NO_EQUIL;
            }
            else
            {
                MSQ_SETERR( err )( "equil flag not set to 0 or 1", MsqError::INVALID_STATE );
            }
        }
    }
    switch( status )
    {
        case MSQ_NO_EQUIL:
            MSQ_PRINT( 3 )( "Not an equilibrium point\n" );
            equil = false;
            break;
        case MSQ_EQUIL:
            MSQ_PRINT( 3 )( "An equilibrium point\n" );
            equil = true;
            break;
        default:
            MSQ_PRINT( 3 )( "Failed to determine equil or not; status = %d\n", status );
    }
    //    FUNCTION_TIMER_END();
    return ( equil );
}

void NonSmoothDescent::form_grammian( const std::vector< Vector3D >& vec, MsqError& /*err*/ )
{
    /* form the grammian with the dot products of the gradients */
    const size_t num_active = mActive.active_ind.size();
    mG.resize( num_active );
    for( size_t i = 0; i < num_active; i++ )
        for( size_t j = i; j < num_active; j++ )
            mG( i, j ) = vec[i] % vec[j];
}

bool NonSmoothDescent::check_equilibrium( OptStatus& status, MsqError& err )
{
    std::vector< Vector3D > dir;

    // TODO - this subroutine is no longer clear to me... is it still
    // appropriate for quads and hexes?  I think it might be in 2D, but
    // 3D is less clear.  Is there a more general algorithm to use?
    // ask Todd/check in numerical optimization

    bool equil              = false;
    const size_t num_active = mActive.active_ind.size();
    ;

    if( num_active == qmHandles.size() )
    {
        equil  = true;<--- equil is assigned
        status = MSQ_EQUILIBRIUM;
        MSQ_PRINT( 3 )( "All the function values are in the active set\n" );
    }

    /* set up the active mGradient directions */
    this->get_active_directions( mGradient, dir, err );
    MSQ_ERRZERO( err );

    /* normalize the active directions */
    for( size_t j = 0; j < num_active; j++ )
        dir[j] /= dir[j].length();

    equil = this->convex_hull_test( dir, err );<--- equil is overwritten
    if( equil ) status = MSQ_EQUILIBRIUM;

    return equil;
}

static double condition3x3( const double A[9] )
{
    //   int ierr;
    double a11, a12, a13;
    double a21, a22, a23;
    double a31, a32, a33;
    //   double s1, s2, s4, s3, t0;
    double s1, s2, s3;
    double denom;
    //   double one = 1.0;
    double temp;
    bool zero_denom = true;

    a11 = A[0];
    a12 = A[1];
    a13 = A[2];
    a21 = A[3];
    a22 = A[4];
    a23 = A[5];
    a31 = A[6];
    a32 = A[7];
    a33 = A[8];

    denom = -a11 * a22 * a33 + a11 * a23 * a32 + a21 * a12 * a33 - a21 * a13 * a32 - a31 * a12 * a23 + a31 * a13 * a22;

    if( ( fabs( a11 ) > EPSILON ) && ( fabs( denom / a11 ) > EPSILON ) )
    {
        zero_denom = false;
    }
    if( ( fabs( a22 ) > EPSILON ) && ( fabs( denom / a22 ) > EPSILON ) )
    {
        zero_denom = false;
    }
    if( ( fabs( a33 ) > EPSILON ) && ( fabs( denom / a33 ) > EPSILON ) )
    {
        zero_denom = false;
    }

    if( zero_denom )
    {
        return HUGE_VAL;
    }
    else
    {
        s1 = sqrt( a11 * a11 + a12 * a12 + a13 * a13 + a21 * a21 + a22 * a22 + a23 * a23 + a31 * a31 + a32 * a32 +
                   a33 * a33 );

        temp = ( -a22 * a33 + a23 * a32 ) / denom;
        s3   = temp * temp;
        temp = ( a12 * a33 - a13 * a32 ) / denom;
        s3 += temp * temp;
        temp = ( a12 * a23 - a13 * a22 ) / denom;
        s3 += temp * temp;
        temp = ( a21 * a33 - a23 * a31 ) / denom;
        s3 += temp * temp;
        temp = ( a11 * a33 - a13 * a31 ) / denom;
        s3 += temp * temp;
        temp = ( a11 * a23 - a13 * a21 ) / denom;
        s3 += temp * temp;
        temp = ( a21 * a32 - a22 * a31 ) / denom;
        s3 += temp * temp;
        temp = ( -a11 * a32 + a12 * a31 ) / denom;
        s3 += temp * temp;
        temp = ( -a11 * a22 + a12 * a21 ) / denom;
        s3 += temp * temp;

        s2 = sqrt( s3 );
        return s1 * s2;
    }
}

double NonSmoothDescent::SymmetricMatrix::condition3x3() const
{
    double values[9] =
        { operator()( 0, 0 ), operator()( 0, 1 ), operator()( 0, 2 ), operator()( 1, 0 ), operator()( 1, 1 ),
          operator()( 1, 2 ), operator()( 2, 0 ), operator()( 2, 1 ), operator()( 2, 2 ) };
    return MBMesquite::condition3x3( values );
}

void NonSmoothDescent::singular_test( int n, const Matrix3D& A, bool& singular, MsqError& err )
{
    //    int test;
    //    double determinant;
    double cond;

    if( ( n > 3 ) || ( n < 1 ) )
    {
        MSQ_SETERR( err )( "Singular test works only for n=1 to n=3", MsqError::INVALID_ARG );
        return;
    }

    singular = true;
    switch( n )
    {
        case 1:
            if( A[0][0] > 0 ) singular = false;
            break;
        case 2:
            if( fabs( A[0][0] * A[1][1] - A[0][1] * A[1][0] ) > EPSILON ) singular = false;
            break;
        case 3:
            /* calculate the condition number */
            cond = condition3x3( A[0] );
            if( cond < 1E14 ) singular = false;
            break;
    }
}

void NonSmoothDescent::form_PD_grammian( MsqError& err )
{
    int i, j, k;
    int g_ind_1;
    bool singular = false;

    const int num_active = mActive.active_ind.size();

    /* this assumes the grammian has been formed */
    for( i = 0; i < num_active; i++ )
    {
        for( j = i; j < num_active; j++ )
        {
            if( mG( i, j ) == -1 )
            {
                MSQ_SETERR( err )( "Grammian not computed properly", MsqError::INVALID_STATE );
                return;
            }
        }
    }

    /* use the first gradient in the active set */
    g_ind_1    = 0;
    mPDG[0][0] = mG( 0, 0 );
    pdgInd[0]  = mActive.active_ind[0];

    /* test the rest and add them as appropriate */
    k = 1;
    i = 1;
    while( ( k < 3 ) && ( i < num_active ) )
    {
        mPDG[0][k] = mPDG[k][0] = mG( 0, i );
        mPDG[k][k]              = mG( i, i );
        if( k == 2 )
        { /* add the dot product of g1 and g2 */
            mPDG[1][k] = mPDG[k][1] = mG( g_ind_1, i );
        }
        this->singular_test( k + 1, mPDG, singular, err );
        if( !singular )
        {
            pdgInd[k] = mActive.active_ind[i];
            if( k == 1 ) g_ind_1 = i;
            k++;
        }
        i++;
    }
}

void NonSmoothDescent::search_edges_faces( const Vector3D* dir, Vector3D& mSearch, MsqError& /*err*/ )
{
    bool viable;
    double a, b, denom;
    Vector3D g_bar;
    Vector3D temp_search( 0, 0, 0 ); /* initialize the search direction to 0,0 */
    double projection, min_projection;

    const size_t num_active = mActive.active_ind.size();
    ;

    /* Check for viable faces */
    min_projection = HUGE_VAL;
    for( size_t i = 0; i < num_active; i++ )
    {
        /* FACE I */
        viable = true;

        /* test the viability */
        for( size_t j = 0; j < num_active; j++ )
        { /* lagrange multipliers>0 */
            if( mG( j, i ) < 0 ) viable = false;
        }

        /* find the minimum of viable directions */
        if( ( viable ) && ( mG( i, i ) < min_projection ) )
        {
            min_projection = mG( i, i );
            temp_search    = dir[i];
            mSteepest      = mActive.active_ind[i];
        }

        /* INTERSECTION IJ */
        for( size_t j = i + 1; j < num_active; j++ )
        {
            viable = true;

            /* find the coefficients of the intersection
               and test the viability */
            denom = 2 * mG( i, j ) - mG( i, i ) - mG( j, j );
            a = b = 0;
            if( fabs( denom ) > EPSILON )
            {
                b = ( mG( i, j ) - mG( i, i ) ) / denom;
                a = 1 - b;
                if( ( b < 0 ) || ( b > 1 ) ) /* 0 < b < 1 */
                    viable = false;
                for( size_t k = 0; k < num_active; k++ )
                { /* lagrange multipliers>0 */
                    if( ( a * mG( k, i ) + b * mG( k, j ) ) <= 0 ) viable = false;
                }
            }
            else
            {
                viable = false; /* Linearly dependent */
            }

            /* find the minimum of viable directions */
            if( viable )
            {
                g_bar      = a * dir[i] + b * dir[j];
                projection = g_bar % g_bar;
                if( projection < min_projection )
                {
                    min_projection = projection;
                    temp_search    = g_bar;
                    mSteepest      = mActive.active_ind[i];
                }
            }
        }
    }
    if( optStatus != MSQ_EQUILIBRIUM )
    {
        mSearch = temp_search;
    }
}

bool NonSmoothDescent::solve2x2( double a11,
                                 double a12,
                                 double a21,
                                 double a22,
                                 double b1,
                                 double b2,
                                 double x[2],
                                 MsqError& /*err*/ )
{
    double factor;

    /* if the system is not singular, solve it */
    if( fabs( a11 * a22 - a21 * a12 ) > EPSILON )
    {
        if( fabs( a11 ) > EPSILON )
        {
            factor = ( a21 / a11 );
            x[1]   = ( b2 - factor * b1 ) / ( a22 - factor * a12 );
            x[0]   = ( b1 - a12 * x[1] ) / a11;
        }
        else if( fabs( a21 ) > EPSILON )
        {
            factor = ( a11 / a21 );
            x[1]   = ( b1 - factor * b2 ) / ( a12 - factor * a22 );
            x[0]   = ( b2 - a22 * x[1] ) / a21;
        }
        return true;
    }
    else
    {
        return false;
    }
}

void NonSmoothDescent::form_reduced_matrix( SymmetricMatrix& P )
{
    const size_t P_size = mActive.active_ind.size() - 1;
    P.resize( P_size );
    for( size_t i = 0; i < P_size; i++ )
    {
        P( i, i ) = mG( 0, 0 ) - 2 * mG( 0, i + 1 ) + mG( i + 1, i + 1 );
        for( size_t j = i + 1; j < P_size; j++ )
        {
            P( i, j ) = mG( 0, 0 ) - mG( 0, j + 1 ) - mG( i + 1, 0 ) + mG( i + 1, j + 1 );
        }
    }
}

void NonSmoothDescent::get_min_estimate( double* final_est, MsqError& /*err*/ )
{
    double est_imp;

    *final_est = -HUGE_VAL;
    for( size_t i = 0; i < mActive.active_ind.size(); i++ )
    {
        est_imp = -mAlpha * mGS[mActive.active_ind[i]];
        if( est_imp > *final_est ) *final_est = est_imp;
    }
    if( *final_est == 0 )
    {
        *final_est = -HUGE_VAL;
        for( size_t i = 0; i < qmHandles.size(); i++ )
        {
            est_imp = -mAlpha * mGS[i];
            if( ( est_imp > *final_est ) && ( fabs( est_imp ) > EPSILON ) )
            {
                *final_est = est_imp;
            }
        }
    }
}

void NonSmoothDescent::get_gradient_projections( const Vector3D& mSearch, MsqError& /*err*/ )
{
    for( size_t i = 0; i < qmHandles.size(); i++ )
        mGS[i] = mGradient[i] % mSearch;

    MSQ_PRINT( 3 )( "steepest in get_gradient_projections %lu\n", (unsigned long)mSteepest );
}

void NonSmoothDescent::compute_alpha( MsqError& /*err*/ )
{
    double steepest_function;
    double steepest_grad;
    double alpha_i;
    double min_positive_value = HUGE_VAL;

    //    FUNCTION_TIMER_START("Compute Alpha");

    MSQ_PRINT( 2 )( "In compute alpha\n" );

    mAlpha = HUGE_VAL;

    steepest_function = mFunction[mSteepest];
    steepest_grad     = mGS[mSteepest];
    for( size_t i = 0; i < qmHandles.size(); i++ )
    {
        /* if it's not active */
        if( i != mSteepest )
        {
            alpha_i = steepest_function - mFunction[i];

            if( fabs( mGS[mSteepest] - mGS[i] ) > 1E-13 )
            {
                /* compute line intersection */
                alpha_i = alpha_i / ( steepest_grad - mGS[i] );
            }
            else
            {
                /* the lines don't intersect - it's not under consideration*/
                alpha_i = 0;
            }
            if( ( alpha_i > minStepSize ) && ( fabs( alpha_i ) < fabs( mAlpha ) ) )
            {
                mAlpha = fabs( alpha_i );
                MSQ_PRINT( 3 )( "Setting alpha %lu %g\n", (unsigned long)i, alpha_i );
            }
            if( ( alpha_i > 0 ) && ( alpha_i < min_positive_value ) )
            {
                min_positive_value = alpha_i;
            }
        }
    }

    if( ( mAlpha == HUGE_VAL ) && ( min_positive_value != HUGE_VAL ) )
    {
        mAlpha = min_positive_value;
    }

    /* if it never gets set, set it to the default */
    if( mAlpha == HUGE_VAL )
    {
        mAlpha = maxAlpha;
        MSQ_PRINT( 3 )( "Setting alpha to the maximum step length\n" );
    }

    MSQ_PRINT( 3 )( "  The initial step size: %f\n", mAlpha );

    //    FUNCTION_TIMER_END();
}

void NonSmoothDescent::print_active_set( const ActiveSet& active_set, const std::vector< double >& func )
{
    if( active_set.active_ind.empty() ) MSQ_DBGOUT( 3 ) << "No active values\n";
    /* print the active set */
    for( size_t i = 0; i < active_set.active_ind.size(); i++ )
    {
        MSQ_PRINT( 3 )
        ( "Active value %lu:   %f \n", (unsigned long)i + 1, func[active_set.active_ind[i]] );
    }
}

void NonSmoothDescent::init_opt( PatchData& pd, MsqError& err )
{
    qmHandles.clear();
    currentQM->get_evaluations( pd, qmHandles, true, err );MSQ_ERRRTN( err );

    MSQ_PRINT( 2 )( "\nInitializing Optimization \n" );

    /* for the purposes of initialization will be set to zero after */
    optStatus = MSQ_NO_STATUS;
    mSteepest = 0;
    mAlpha    = 0;
    maxAlpha  = 0;

    MSQ_PRINT( 3 )( "  Initialized Constants \n" );
    pdgInd[0] = pdgInd[1] = pdgInd[2] = -1;
    mPDG                              = Matrix3D( 0.0 );

    MSQ_PRINT( 3 )( "  Initialized search and PDG \n" );
    mFunction.clear();
    mFunction.resize( qmHandles.size(), 0.0 );
    testFunction.clear();
    testFunction.resize( qmHandles.size(), 0.0 );
    originalFunction.clear();
    originalFunction.resize( qmHandles.size(), 0.0 );
    mGradient.clear();
    mGradient.resize( qmHandles.size(), Vector3D( 0.0 ) );
    mGS.resize( qmHandles.size() );

    MSQ_PRINT( 3 )( "  Initialized function/gradient \n" );
    mG.resize( qmHandles.size() );
    mG.fill( -1 );
    MSQ_PRINT( 3 )( "  Initialized G\n" );

    prevActiveValues.clear();
    prevActiveValues.reserve( 32 );
    MSQ_PRINT( 3 )( "  Initialized prevActiveValues\n" );
}

void NonSmoothDescent::init_max_step_length( PatchData& pd, MsqError& err )
{
    size_t i, j;
    double max_diff = 0;
    double diff     = 0;<--- Variable 'diff' is assigned a value that is never used.

    MSQ_PRINT( 2 )( "In init_max_step_length\n" );

    /* check that the input data is correct */
    if( pd.num_elements() == 0 )
    {
        MSQ_SETERR( err )( "Num incident vtx = 0\n", MsqError::INVALID_MESH );
        return;
    }

    /* find the maximum distance between two incident vertex locations */
    const MsqVertex* coords = pd.get_vertex_array( err );
    for( i = 0; i < pd.num_nodes() - 1; i++ )
    {
        for( j = i; j < pd.num_nodes(); j++ )
        {
            diff = ( coords[i] - coords[j] ).length_squared();
            if( max_diff < diff ) max_diff = diff;
        }
    }
    max_diff = sqrt( max_diff );
    if( max_diff == 0 )
    {
        MSQ_SETERR( err )
        ( "Maximum distance between incident vertices = 0\n", MsqError::INVALID_MESH );
        return;
    }
    maxAlpha = max_diff / 100;

    MSQ_PRINT( 3 )( "  Maximum step is %g\n", maxAlpha );
}

}  // namespace MBMesquite