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
//! Implementation of the HEALPix framework.  
//! See papers:  
//! * Gorsky2005: "HEALPix: A Framework for High-Resolution Discretization and Fast Analysis of Data
//!                Distributed on the Sphere", Górski, K. M. et al., 2005; 2005ApJ...622..759G.
//! * Calabretta2004: "Mapping on the HEALPix grid", Calabretta, M. R., 2004; 2004astro.ph.12607C
//! * Calabretta2007: "Mapping on the HEALPix grid", Calabretta, M. R. et Roukema, B. F., 2007; 2007MNRAS.381..865C
//! * Reinecke2015: "Efficient data structures for masks on 2D grids",  Reinecke, M. et Hivon, E., 2015; 2015A&A...580A.132R

#![cfg_attr(test, feature(test))]

#[cfg(test)]
extern crate test;

use std::sync::Once;
use std::f64::consts::{PI};

const SQRT2: f64 = 1.41421356237309504880_f64;
const SQRT6: f64 = 2.44948974278317809819_f64;
const ONE_OVER_SQRT6: f64 = 0.40824829046386301636_f64;
const HALF: f64 = 0.5_f64;

/// Upper limit on sqrt(3(1-|z|)) to consider that we are not near from the poles
const EPS_POLE: f64 = 1e-13_f64;

/// Constant = pi/2.
/// 
/// ```rust
/// use cdshealpix::{HALF_PI};
/// use std::f64::consts::PI;
/// assert_eq!(PI / 2f64, HALF_PI);
/// ```
pub const HALF_PI: f64 = HALF * PI;

/// Constant = 4/pi.
/// 
/// ```rust
/// use cdshealpix::{FOUR_OVER_PI};
/// use std::f64::consts::PI;
/// assert_eq!(4f64 / PI, FOUR_OVER_PI);
/// ```
pub const FOUR_OVER_PI: f64 = 4_f64 / PI;

/// Constant = pi/4.
/// 
/// ```rust
/// use cdshealpix::{PI_OVER_FOUR};
/// use std::f64::consts::PI;
/// assert_eq!(PI / 4f64, PI_OVER_FOUR);
/// ```
pub const PI_OVER_FOUR: f64 = 0.25_f64 * PI;

/// Constant = 29, i.e. the largest possible depth we can store on a signed positive long
/// (4 bits for base cells + 2 bits per depth + 2 remaining bits (1 use in the unique notation).
/// 
/// ```rust
/// use cdshealpix::{DEPTH_MAX};
/// assert_eq!(29, DEPTH_MAX);
/// ```
pub const DEPTH_MAX: u8 = 29;

/// Constant = nside(29), i.e. the largest possible nside available when we store HEALPix hash
/// on a u64.
///
/// ```rust
/// use cdshealpix::{DEPTH_MAX, NSIDE_MAX, nside};
/// assert_eq!(nside(DEPTH_MAX), NSIDE_MAX);
/// ```
pub const NSIDE_MAX: u32 = 536870912;

/// Limit on the latitude (in radians) between the equatorial region and the polar caps.
/// Equals asin(2/3) = 0.7297276562269663 radians ~= 41,81 degrees.
/// Written $\theta_X$ in Calabretta2007.
/// 
/// ```rust
/// use cdshealpix::{TRANSITION_LATITUDE};
/// assert_eq!(f64::asin(2f64 / 3f64), TRANSITION_LATITUDE);
/// ```
pub const TRANSITION_LATITUDE: f64 = 0.72972765622696636344_f64; // asin(2/3)
/// Limit on |z|=|sin(lat)| between the equatorial region and the polar caps.
/// Equals 2/3, see Eq. (1) in Gorsky2005.
pub const TRANSITION_Z: f64 = 2_f64 / 3_f64;
const ONE_OVER_TRANSITION_Z: f64 = 1.5_f64;

const F64_SIGN_BIT_MASK: u64 = 0x8000000000000000;
const F64_BUT_SIGN_BIT_MASK: u64 = 0x7FFFFFFFFFFFFFFF;

/// For each HEALPix depth, stores the smallest distance from an edge of a cell to the opposite
/// edge of the same cell. If the radius of a cone is smaller than this distance, we know that
/// it will overlap maximum 9 pixels (the pixel containing the center of the cone plus
/// the 8 neighbours).  
/// In practice, this distance if the distance between the point of coordinate
/// (0, TRANSITION_LATITUDE) and it nearest point on the Northeast edge of the
/// cell of base hash 0 and coordinates in the base hash (x=0, y=nside-1).  
/// IMPORTANT REMARK:  
/// - this value is larger than the smallest center to vertex distance
/// - this value x2 is larger than the smallest diagonal (NS or EW)
/// - this value x2 is larger than the smallest edge
/// - BUT there is no case in which the value is larger than the four center-to-vertex distance
/// - BUT there is no case in which the value x2 is larger than both diagonals
/// => this radius is smaller than the smaller circumcircle radius (=> no cone having the smaller
/// -edge-to-opposite-edge-radius radius can contains the 4 vertices of a cell (but 3 is ok)
/// vertices
static SMALLER_EDGE2OPEDGE_DIST: [f64; 30] =  [
    0.8410686705685088,    // depth = 0
    0.37723631722170053,   // depth = 1
    0.18256386461918295,   // depth = 2
    0.09000432499034523,   // depth = 3
    0.04470553761855741,   // depth = 4
    0.02228115704023076,   // depth = 5
    0.011122977211214961,  // depth = 6
    0.005557125022105058,  // depth = 7
    0.0027774761500209185, // depth = 8
    0.0013884670480328143, // depth = 9
    6.941658374603201E-4,  // depth = 10
    3.4706600585087755E-4, // depth = 11
    1.7352877579970442E-4, // depth = 12
    8.676333125510362E-5,  // depth = 13
    4.338140148342286E-5,  // depth = 14
    2.1690634707822447E-5, // depth = 15
    1.084530084565172E-5,  // depth = 16
    5.422646295795749E-6,  // depth = 17
    2.711322116099695E-6,  // depth = 18
    1.3556608000873442E-6, // depth = 19
    6.778303355805395E-7,  // depth = 20
    3.389151516386149E-7,  // depth = 21
    1.69457571754776E-7,   // depth = 22
    8.472878485272006E-8,  // depth = 23
    4.236439215502565E-8,  // depth = 24
    2.1182195982014308E-8, // depth = 25
    1.0591097960375205E-8, // depth = 26
    5.295548939447981E-9,  // depth = 27
    2.647774429917369E-9,  // depth = 28
    1.3238871881399636E-9  // depth = 29
];

// Idee pour cone: surface du cone => ncells !! (facteur variant 8, 4, 2, ... 1.2 

/// Latitude, in the equatorial region, for which the distance from the cell center to its four
/// vertices is almost equal on the sky (i.e. the shape of the cell on the sky is close to a square).
/// The larger the depth, the better the approximation (based on differential calculus).
/// > dX = dY = 1 / nside (center to vertex distance)
/// > X = 4/pi * lon     => dX = 4/pi dlon
/// > Y = 3/2 * sin(lat) => dY = 3/2 * cos(lat) dlat
/// > dlon * cos(lat) = dlat (same distance on the sky)
/// > => cos^2(lat) = 2/3 * 4/pi
/// > => lat = arccos(sqrt(2/3 * 4/pi)) ~= 22.88 deg ~= 0.39934 rad
/// 
/// ```rust
/// use cdshealpix::{TRANSITION_Z, FOUR_OVER_PI, LAT_OF_SQUARE_CELL};
/// assert!(f64::abs(f64::acos(f64::sqrt(TRANSITION_Z * FOUR_OVER_PI)) - LAT_OF_SQUARE_CELL) < 1e-15_f64);
/// ```
pub static LAT_OF_SQUARE_CELL: f64 = 0.39934019947897773410_f64;
/// Simply the consine of LAT_OF_SQUARE_CELL
static COS_LAT_OF_SQUARE_CELL: f64 = 0.92131773192356127804_f64;

/// Array storing pre-computed values for each of the 30 possible depth (from 0 to 29)
/// Info: I would have prefered to compute those quantities at compilation time, and thus have
/// a `static CSTS_C2V: [ConstantsC2V; 30]`. Unfortunately:
/// - macro do not seems work with static arrays
/// - const fn is not stable and can only use const fn (so no min/max/sin/ ...) :o/
static mut CSTS_C2V: [Option<ConstantsC2V>; 30] = [
  None, None, None, None, None, None, None, None, None, None, None, None, None, None, None,
  None, None, None, None, None, None, None, None, None, None, None, None, None, None, None
];
// Found here: https://stackoverflow.com/questions/28656387/initialize-a-large-fixed-size-array-with-non-copy-types
// I wanted to use it to set a static array. So far it is not possible.
// I hope in the future to be able to compute such an array at compilation time.
/*macro_rules! make_array {
  ($n: expr, $constructor: expr) => {
    {
      let mut items: [_; $n] = mem::uninitialized();
      for (i, place) in items.iter_mut().enumerate() {
          ptr::write(place, $constructor(i as u8));
      }
      items
    }
  }
}
static CSTS_C2V: [ConstantsC2V; 30] = unsafe { make_array!(30, |depth| new_cst_c2v(depth)) };
*/

/*
/// Defines a simple range `from` a given index (inclusive) `to` a given index (exclusive) 
pub struct Range {
  pub from: u64,
  pub to: u64,
}

impl Range {
  
  pub fn new(from: u64, to: u64) -> Range {
    if from >= to {
      panic!("Illegal argument: from '{}' must be lower than to '{}'.", from, to);
    }
    new_unsafe(from, to)
  }
  
  fn new_unsafe(from: u64, to: u64) -> Range {
    Range{from, to}
  }
  
}
*/

/// See the get_or_create function, each object is used for the lazy instantiation of the 
/// layer of the corresponding depht.
/// Info: Unfortunatly Default::default(); do no work with static arrays :o/
static CSTS_C2V_INIT: [Once; 30] = [
  Once::new(), Once::new(), Once::new(), Once::new(), Once::new(), Once::new(), Once::new(),
  Once::new(), Once::new(), Once::new(), Once::new(), Once::new(), Once::new(), Once::new(),
  Once::new(), Once::new(), Once::new(), Once::new(), Once::new(), Once::new(), Once::new(),
  Once::new(), Once::new(), Once::new(), Once::new(), Once::new(), Once::new(), Once::new(),
  Once::new(), Once::new()
];
/// Lazy factory method: instantiate a new Layer at the first call for a given depth; after the 
/// first call, returns an already instantiated Layer.
/// # Info
/// This method resort to a double-checked lock, ensuring thread-safety.
fn get_or_create(depth: u8) -> &'static ConstantsC2V {
  unsafe {
    // Inspired from the Option get_or_insert_with method, modified to ensure thread safety with
    // https://doc.rust-lang.org/std/sync/struct.Once.html
    // This implements a double-checked lock
    match CSTS_C2V[depth as usize] {
      Some(ref v) =>  return v,
      None => {
        CSTS_C2V_INIT[depth as usize].call_once(|| {
          CSTS_C2V[depth as usize] = Some(ConstantsC2V::new(depth));
        });
      },
    }
    match CSTS_C2V[depth as usize] {
      Some(ref v) => v,
      _ => unreachable!(),
    }
  }
}

struct ConstantsC2V {
  slope_npc: f64,
  intercept_npc: f64,
  slope_eqr: f64,
  intercept_eqr: f64,
  coeff_x2_eqr: f64,
  coeff_cst_eqr: f64,
}

impl ConstantsC2V {
 fn new(depth: u8) -> ConstantsC2V {
    let nside = nside_unsafe(depth);
    let dist_cw = 1.0_f64 / (nside as f64); // Center to West (or East) vertex distance on the transition latitude
    let one_min_dist_cw = 1.0_f64 - dist_cw;
    // NPC, see comment of function largest_c2v_dist_in_npc()
    let lat_north = f64::asin(1_f64 - (pow2(one_min_dist_cw) / 3_f64));
    let mut d_min = lat_north - TRANSITION_LATITUDE;
    let mut d_max: f64 = sphe_dist(
      squared_half_segment(
        PI_OVER_FOUR * dist_cw, d_min,
        lat_north.cos(),  TRANSITION_LATITUDE.cos()));
    // - linear approx
    let slope_npc: f64 = (d_max - d_min) / (PI_OVER_FOUR * one_min_dist_cw); // a = (yB - yA) / (xB - xA); with xA = 0
    let intercept_npc: f64 = d_min;                                          // b = yA - a * xA          ; with xA = 0
    // EQR TOP, see comment of function largest_c2v_dist_in_eqr_top()
    d_min =  FOUR_OVER_PI * dist_cw * COS_LAT_OF_SQUARE_CELL;
    d_max = TRANSITION_LATITUDE - f64::asin(one_min_dist_cw * TRANSITION_Z);
    let slope_eqr: f64 = (d_max - d_min) / (TRANSITION_LATITUDE - LAT_OF_SQUARE_CELL); // a = (yB - yA) / (xB - xA) 
    let intercept_eqr: f64 = d_min - slope_eqr * LAT_OF_SQUARE_CELL;                   // b = yA - a * xA 
    // EQR BOTTOM,  see comment of function largest_c2v_dist_in_eqr_bottom()
    d_max = FOUR_OVER_PI * dist_cw;
    let coeff_cst_eqr: f64 = d_max;
    let coeff_x2_eqr: f64 = (d_min - d_max) / pow2(LAT_OF_SQUARE_CELL);
    // Struct creation
    ConstantsC2V {
      slope_npc,
      intercept_npc,
      slope_eqr,
      intercept_eqr,
      coeff_x2_eqr,
      coeff_cst_eqr
    }
  }
}

#[inline]
fn haversine_dist(p1_lon: f64, p1_lat: f64, p2_lon: f64, p2_lat: f64) -> f64 {
  let shs = squared_half_segment(
    p2_lon - p1_lon, p2_lat - p1_lat, 
    p1_lat.cos(), p2_lat.cos());
  sphe_dist(shs)
}

/// Returns the angular distance corresponding to the given squared half great-circle arc segment
#[inline]
fn sphe_dist(squared_half_segment: f64) -> f64 {
  squared_half_segment.sqrt().asin().twice()
}

/// Returns `(s/2)^2` with `s` the segment (i.e. the Euclidean distance) between 
/// the two given points  `P1` and `P2` on the unit-sphere.
/// We recall that `s = 2 sin(ad/2)` with `ad` the angular distance between the two points.
/// # Input
/// - `dlon` the longitude difference, i.e. (P2.lon - P1.lon), in radians
/// - `dlat` the latitude difference, i.e. (P2.lat - P1.lat), in radians
/// - `cos_lat1` cosine of the latitude of the first point
/// - `cos_lat2` cosine of the latitude of the second point
#[inline]
fn squared_half_segment(dlon: f64, dlat: f64, cos_lat1: f64, cos_lat2: f64) -> f64 {
  dlat.half().sin().pow2() + cos_lat1 * cos_lat2 * dlon.half().sin().pow2()
}


#[inline]
fn to_squared_half_segment(spherical_distance: f64) -> f64 {
  spherical_distance.half().sin().pow2()
}

/*#[inline]
fn to_squared_half_segments(spherical_distances: &mut [f64]) {
  spherical_distances.iter_mut().for_each(|d| *d = to_squared_half_segment(*d));
}*/

#[inline]
fn pow2(x: f64) -> f64 {
  x * x
}

/*impl f64 {
  #[inline]
  pub fn pow2(self) -> f64 {
    self * self
  }
}*/

/// Simple trait used to implements `pow2`, `twice` and `half` on f64.
pub trait Customf64 {
  #[inline]
  fn pow2(self) -> f64;
  #[inline]
  fn twice(self) -> f64;
  #[inline]
  fn half(self) -> f64;
  #[inline]
  fn div_eucl(self, rhs: f64) -> f64;
}

impl Customf64 for f64 {
  /// Returns x^2
  #[inline]
  fn pow2(self) -> f64 {
    self * self // or powi ?
  }
  /// Returns 2 * x
  #[inline]
  fn twice(self) -> f64 {
    2.0 * self // self + self (I hope the compiler know the simple shift bit to be used for x2)
  }
  /// Returns x / 2
  #[inline]
  fn half(self) -> f64 {
    0.5 * self
  }

  /// [Duplicated code](https://doc.rust-lang.org/std/primitive.f64.html#method.div_euc), because
  /// it is unstable so far.
  #[inline]
  fn div_eucl(self, rhs: f64) -> f64 {
    let q = (self / rhs).trunc();
    if self % rhs < 0.0 {
      return if rhs > 0.0 { q - 1.0 } else { q + 1.0 }
    }
    q
  }
}

/// All types that implement `Write` get methods defined in `WriteBytesExt`
/// for free.
// impl<F: f64> Customf64 for F {}*/

/// Returns an upper limit on the distance between a cell center around the given position 
/// and its furthest vertex.
/// # Params
/// - `depth` the depth of the cell
/// - `lon` the longitude of the point on the unit-sphere, in radians 
/// - `lat` the latitude of the point on the unit-sphere, in radians
/// 
/// # Result
/// The following plot shows, for the depth 8, the real largest distances (in red) and the result
/// of this method (in blue).  
/// WARNING: the units of `(lon, lat)` on the plot are *degrees*, while the distance is in *mas*  
/// Credit: plot made using [TOPCAT](http://www.star.bris.ac.uk/~mbt/topcat/)
/// ![CenterToVertexDist](https://raw.githubusercontent.com/cds-astro/cds-healpix-rust/master/resources/4doc/d_center_vertex.png)
/// 
pub fn largest_center_to_vertex_distance(depth: u8, lon: f64, lat: f64) -> f64 {
  // Specific case for depth 0
  if depth == 0 {
    return HALF_PI - TRANSITION_LATITUDE;
  }
  // Regular case
  let lat_abs = lat.abs();
  if lat_abs >= TRANSITION_LATITUDE {
    largest_c2v_dist_in_npc(lon, get_or_create(depth))
  } else if lat_abs >= LAT_OF_SQUARE_CELL {
    largest_c2v_dist_in_eqr_top(lat_abs, get_or_create(depth))
  } else {
    largest_c2v_dist_in_eqr_bottom(lat_abs, get_or_create(depth))
  }
}

/// Returns an upper limit on the distance between a cell center and it furthest vertex, for
/// all the cells in the region covered by a cone of given center and radius.  
/// It is an extension of [largest_center_to_vertex_distance](#fn.largest_center_to_vertex_distance)
/// # Params
/// - `depth` the depth of the cell
/// - `lon` the longitude of the point on the unit-sphere, in radians 
/// - `lat` the latitude of the point on the unit-sphere, in radians
/// - `radius` the radius of the cone, in radians
///
pub fn largest_center_to_vertex_distance_with_radius(depth: u8, lon: f64, lat: f64, radius: f64) -> f64 {
  // Specific case for depth 0
  if depth == 0 {
    return HALF_PI - TRANSITION_LATITUDE;
  }
  // Regular case
  let lat_abs = lat.abs();
  let lat_max = lat_abs + radius;
  let lat_min = lat_abs - radius;
  if lat_max >= TRANSITION_LATITUDE {
    largest_c2v_dist_in_npc_with_radius(lon, radius, get_or_create(depth))
  } else if lat_min >= LAT_OF_SQUARE_CELL {
    largest_c2v_dist_in_eqr_top_with_radius(lat_abs, radius, get_or_create(depth))
  } else if lat_max <= LAT_OF_SQUARE_CELL {
    largest_c2v_dist_in_eqr_bottom_with_radius(lat_abs, radius, get_or_create(depth))
  } else {
    let csts = get_or_create(depth);
    f64::max(
      largest_c2v_dist_in_eqr_top_with_radius(lat_abs, radius, csts),
      largest_c2v_dist_in_eqr_bottom_with_radius(lat_abs, radius, csts)
    )
  }
}

/// Same as [largest_center_to_vertex_distance_with_radius](#fn.largest_center_to_vertex_distance_with_radius)
/// but making the computation for several depths at the same time.
pub fn largest_center_to_vertex_distances_with_radius(mut from_depth: u8, to_depth: u8, lon: f64, lat: f64, radius: f64) -> Box<[f64]> {
  let mut vec: Vec<f64> = Vec::with_capacity((to_depth - from_depth) as usize);
  // Specific case for depth 0
  if from_depth == 0 {
    vec.push(HALF_PI - TRANSITION_LATITUDE);
    from_depth = 1_u8;
  }
  // Regular case
  let lat_abs = lat.abs();
  let lat_max = lat_abs + radius;
  let lat_min = lat_abs - radius;
  if lat_max >= TRANSITION_LATITUDE {
    let mut lon = (PI_OVER_FOUR - (lon % HALF_PI)).abs();
    lon = f64::min(lon + radius, PI_OVER_FOUR);
    for depth in from_depth..to_depth {
      let csts = get_or_create(depth);
      vec.push(linear_approx(lon, csts.slope_npc, csts.intercept_npc));
    }
  } else if lat_min >= LAT_OF_SQUARE_CELL {
    for depth in from_depth..to_depth {
      vec.push(largest_c2v_dist_in_eqr_top(lat_max, get_or_create(depth)));
    }
  } else if lat_max <= LAT_OF_SQUARE_CELL {
    let val_min = f64::max(lat_min, 0_f64);
    for depth in from_depth..to_depth {
      vec.push(largest_c2v_dist_in_eqr_bottom(val_min, get_or_create(depth))); 
    }
  } else {
    let val_max = f64::min(lat_max, TRANSITION_LATITUDE);
    let val_min = f64::max(lat_min, 0_f64);
    for depth in from_depth..to_depth {
      let csts = get_or_create(depth);
      vec.push(
        f64::max(
          largest_c2v_dist_in_eqr_top(val_max, csts),
          largest_c2v_dist_in_eqr_bottom(val_min, csts)
        )
      );
    }
  }
  vec.into_boxed_slice()
}

/// Returns an upper limit on distance between the center of a cell and its furthest vertex.  
/// We assumes that the cell center is located in the North polar cap region (OR IS ON THE
/// TRANSITION LATITUDE).  
/// We use a linear upper limit based on the longitude.
/// - At the transition latitude, we note CN the distance at a base cell border between:
///   - Cell center (C): (x_c = 1/nside, y_c = 1) => (lon_c = pi/4 * 1/nside, lat_c = TRANSITION_LATITUDE)
///   - North vertex (N): (x_n = 0, y_c = 2 - (1 + 1/nside)) => (lon_n = 0, lat_n = asin(1 - ((1 - 1/nside)^2 / 3)) 
///   - using the Haversine formula and SC defined below:
/// > dMax = CN = 2 * asin(sqrt( sin^2(SC/2)  + sin^2(pi/8) * cos(pi/4 * 1/nside)))
/// - At the transition latitude, we not SC the distance at a base cell center between:
///   - South vertex (S): (lon_s = pi/4, lat_s = TRANSITION_LATITUDE)
///   - Cell center: (lon_c = pi/4,  lat_c = asin(1 - ((1 - 1/nside)^2 / 3))
/// > dMin = SC =  asin(1 - ((1 - 1/nside)^2 / 3)) - TRANSITION_LATITUDE
/// - finally, linear approx:
/// > d = ((lon % pi/2) - pi/4) * (dMax - dMin)/(pi/4 * (1 - 1/nside)) + dMin
#[inline]
fn largest_c2v_dist_in_npc(lon: f64, csts: &ConstantsC2V) -> f64 {
  let lon = (PI_OVER_FOUR - (lon % HALF_PI)).abs();
  debug_assert!(0_f64 <= lon && lon <= PI_OVER_FOUR);
  linear_approx(lon, csts.slope_npc, csts.intercept_npc)
}
/// Same as the above method, but taking into account an additional radius
#[inline]
fn largest_c2v_dist_in_npc_with_radius(lon: f64, radius: f64, csts: &ConstantsC2V) -> f64 {
  debug_assert!(0_f64 < radius);
  let mut lon = (PI_OVER_FOUR - (lon % HALF_PI)).abs();
  debug_assert!(0_f64 <= lon && lon <= PI_OVER_FOUR);
  lon = f64::min(lon + radius, PI_OVER_FOUR);
  linear_approx(lon, csts.slope_npc, csts.intercept_npc)
}

/// Returns an upper limit on distance between the center of a cell and its furthest vertex.  
/// We assumes that the cell center is located in the equatorial region,
/// above the latitude at which cells are squares.
/// We use a linear upper limit based on the latitude.
/// - At the latitude in which cells are squares (lat = LAT_OF_SQUARE_CELL), we note CE the 
///   Center-to-East distance: CE (=CW) =  4/pi * 1/nside * cos(LAT_OF_SQUARE_CELL)
/// > dMin = 4/pi * 1/nside * cos(LAT_OF_SQUARE_CELL)
/// - At the latitude under the transition latitude, we note CN the distance with N on the 
///   transition latitude and C at latitude: y_c = 1 - 1/nside => lat_c = arcsin(y_c * 2/3)
/// > dMax = TRANSITION_LATITUDE - arcsin(y_c * 2/3)
/// - finally, linear approx:
/// > d = (lat - LAT_OF_SQUARE_CELL) * (dMax - dMin)/(TRANSITION_LATITUDE - LAT_OF_SQUARE_CELL) + dMin
#[inline]
fn largest_c2v_dist_in_eqr_top(lat_abs: f64, csts: &ConstantsC2V) -> f64 {
  debug_assert!(LAT_OF_SQUARE_CELL <= lat_abs && lat_abs < TRANSITION_LATITUDE);
  linear_approx(lat_abs, csts.slope_eqr, csts.intercept_eqr)
}
/// Same as the above method, but taking into account an additional radius
#[inline]
fn largest_c2v_dist_in_eqr_top_with_radius(lat_abs: f64, radius: f64, csts: &ConstantsC2V) -> f64 {
  debug_assert!(0_f64 < radius);
  debug_assert!(LAT_OF_SQUARE_CELL <= lat_abs && lat_abs < TRANSITION_LATITUDE);
  largest_c2v_dist_in_eqr_top(f64::min(lat_abs + radius, TRANSITION_LATITUDE), csts)
}

/// Returns an upper limit on distance between the center of a cell and its furthest vertex.  
/// We assumes that the cell center is located in the equatorial region,
/// bellow the latitude at which cells are squares.
/// We use a parabola approximation upper limit based on the latitude.
/// - At lat = 0, d_max = pi/4 * 1/nside
/// - At lat = LAT_OF_SQUARE_CELL, d_min = 4/pi * 1/nside * cos(LAT_OF_SQUARE_CELL)
/// - Parabola approx: a * lat^2 + b = dist
///   - At lat = 0, b = d_max
///   - At lat = LAT_OF_SQUARE_CELL, a * LAT_OF_SQUARE_CELL^2 + d_max = d_min
///   - => a = (d_min - d_max) / LAT_OF_SQUARE_CELL^2
///   - => a = d_max (cos(LAT_OF_SQUARE_CELL) - 1) / LAT_OF_SQUARE_CELL^2
///   - => dist =  d_max * (1 + (cos(LAT_OF_SQUARE_CELL) - 1) / LAT_OF_SQUARE_CELL^2 * lat^2)
///   - => dist =  d_max * (1 - (1 - cos(LAT_OF_SQUARE_CELL)) / LAT_OF_SQUARE_CELL^2 * lat^2)
#[inline]
fn largest_c2v_dist_in_eqr_bottom(lat_abs: f64, csts: &ConstantsC2V) -> f64 {
  debug_assert!(0_f64 <= lat_abs && lat_abs <= LAT_OF_SQUARE_CELL);
  csts.coeff_x2_eqr * pow2(lat_abs) + csts.coeff_cst_eqr
}
/// Same as the above method, but taking into account an additional radius.
#[inline]
fn largest_c2v_dist_in_eqr_bottom_with_radius(lat_abs: f64, radius: f64, csts: &ConstantsC2V) -> f64 {
  debug_assert!(0_f64 < radius);
  debug_assert!(0_f64 <= lat_abs && lat_abs <= LAT_OF_SQUARE_CELL);
  largest_c2v_dist_in_eqr_bottom(f64::max(lat_abs - radius, 0_f64), csts)
}

#[inline]
fn linear_approx(x: f64, slope: f64, intercept: f64) -> f64 {
  slope * x + intercept
}

/// Returns, for the given depth, the number of cells along both axis of a base-resolution cell.
///
/// # Input
/// - `depth` must be in `[0, 29]`
///
/// # Output
/// - `nside` = 2^`depth`
///
/// # Panics
/// If `depth` is not valid (see [is_depth](fn.is_depth.html)), this method panics.
///
/// # Examples
///
/// ```rust
/// use cdshealpix::{nside};
///
/// assert_eq!(1, nside(0));
/// assert_eq!(2, nside(1));
/// assert_eq!(4, nside(2));
/// assert_eq!(8, nside(3));
/// assert_eq!(16, nside(4));
/// assert_eq!(32, nside(5));
/// assert_eq!(64, nside(6));
/// assert_eq!(128, nside(7));
/// assert_eq!(256, nside(8));
/// assert_eq!(512, nside(9));
/// assert_eq!(1024, nside(10));
/// assert_eq!(2048, nside(11));
/// assert_eq!(4096, nside(12));
/// assert_eq!(8192, nside(13));
/// assert_eq!(16384, nside(14));
/// assert_eq!(32768, nside(15));
/// assert_eq!(65536, nside(16));
/// assert_eq!(131072, nside(17));
/// assert_eq!(262144, nside(18));
/// assert_eq!(524288, nside(19));
/// assert_eq!(1048576, nside(20));
/// assert_eq!(2097152, nside(21));
/// assert_eq!(4194304, nside(22));
/// assert_eq!(8388608, nside(23));
/// assert_eq!(16777216, nside(24));
/// assert_eq!(33554432, nside(25));
/// assert_eq!(67108864, nside(26));
/// assert_eq!(134217728, nside(27));
/// assert_eq!(268435456, nside(28));
/// assert_eq!(536870912, nside(29));
/// // Using a for loop...
/// for depth in 0..29 {
///     assert_eq!(2u32.pow(depth), nside(depth as u8));
/// }
/// ```
#[inline]
pub fn nside(depth: u8) -> u32 {
    check_depth(depth);
    nside_unsafe(depth)
}

/// Same as [nside](fn.nside.html) except that this version does not check the argument, and thus
/// does not panics if the argument is illegal.
#[inline]
pub fn nside_unsafe(depth: u8) -> u32 {
    1_u32 << depth
}


/// Returns, for the given difference of depth, the number of cells small cells the large cell
/// contains. If the small cell level is 0, the result is the sqaured nside.
///
/// # Input
/// - `delta_depth` must be in `[0, 29]`
///
/// # Output
/// - `nside^2` = 2^2*`delta_depth`
///
/// # Panics
/// If `delta_depth` is not in `[0, 29]`.
///
/// # Examples
///
/// ```rust
/// use cdshealpix::{nside_square};
/// 
/// for delta_depth in 0..29_u8 {
///     assert_eq!(2u64.pow(2 * delta_depth as u32), nside_square(delta_depth));
/// }
/// ```
#[inline]
pub fn nside_square(delta_depth: u8) -> u64 {
  check_depth(delta_depth);
  nside_square_unsafe(delta_depth)
}

/// Same as [nside_square](fn.nside_square.html) except that this version does not check the argument, 
/// and thus does not panics if the argument is illegal.
#[inline]
pub fn nside_square_unsafe(delta_depth: u8) -> u64 {
  1_u64 << (delta_depth << 1)
}

#[inline]
fn check_depth(depth: u8) {
    assert!(is_depth(depth), "Expected depth in [0, 29]");
}

/// Returns `true` if the given argument is a valid depth, i.e. if it is <= [DEPTH_MAX](constant.DEPTH_MAX.html). 
#[inline]
pub fn is_depth(depth: u8) -> bool {
    depth <= DEPTH_MAX
}

/// Returns, for the given `nside`, the number of subdivision of a base-resolution cell (i.e. the depth).
/// 
/// # Input
/// - `nside` must be a power of 2 in `[0, 2^29]`
///
/// # Output
/// - `depth` = `log2(nside)`
///
/// # Panics
/// If `nside` is not valid (see [is_nside](fn.is_nside.html)), this method panics.
///
/// # Examples
///
/// ```rust
/// use cdshealpix::{nside, depth};
///
/// for d in 0..29 {
///     assert_eq!(d, depth(nside(d as u8)));
/// }
/// ```
#[inline]
pub fn depth(nside: u32) -> u8 {
    check_nside(nside);
    depth_unsafe(nside)
}

/// Same as [depth](fn.depth.html) except that this version does not check the argument, and thus
/// does not panics if the argument is illegal.
#[inline]
pub fn depth_unsafe(nside: u32) -> u8 {
    nside.trailing_zeros() as u8
}

#[inline]
fn check_nside(nside: u32) {
    assert!(is_nside(nside), "Nside must be a power of 2 in [1-2^29]");
}

/// Returns `true` if the given argument is a valid `nside`, i.e. 
/// if it is a power of 2, is != 0 and is <= [NSIDE_MAX](constant.NSIDE_MAX.html). 
#[inline]
pub fn is_nside(nside: u32) -> bool {
    is_pow_of_2(nside) && nside > 0 && nside <= NSIDE_MAX
}

/// Determines if an integer is a power of two, including 0.
/// Taken from the "Bit Twiddling Hacks" web page of Sean Eron Anderson.
#[inline]
fn is_pow_of_2(x: u32) -> bool {
    (x & (x - 1)) == 0
}

/// Returns the number of distinct hash value (the number of cells or pixel the unit sphere is
/// devided in) at the given `depth`.
/// 
/// # Input
/// - `depth` must be in `[0, 29]`
///
/// # Output
/// - `n_hash` = `12 * nside^2`
///
/// # Panics
/// If `depth` is not valid (see [is_depth](fn.is_depth.html)), this method panics.
///
/// # Examples
/// 
/// ```rust
/// use cdshealpix::{n_hash};
/// 
/// assert_eq!(12u64, n_hash(0u8));
/// assert_eq!(48u64, n_hash(1u8));
/// assert_eq!(192u64, n_hash(2u8));
/// assert_eq!(768u64, n_hash(3u8));
/// assert_eq!(3072u64, n_hash(4u8));
/// assert_eq!(12288u64, n_hash(5u8));
/// assert_eq!(49152u64, n_hash(6u8));
/// assert_eq!(196608u64, n_hash(7u8));
/// assert_eq!(786432u64, n_hash(8u8));
/// assert_eq!(3145728u64, n_hash(9u8));
/// assert_eq!(12582912u64, n_hash(10u8));
/// assert_eq!(50331648u64, n_hash(11u8));
/// assert_eq!(201326592u64, n_hash(12u8));
/// assert_eq!(805306368u64, n_hash(13u8));
/// assert_eq!(3221225472u64, n_hash(14u8));
/// assert_eq!(12884901888u64, n_hash(15u8));
/// assert_eq!(51539607552u64, n_hash(16u8));
/// assert_eq!(206158430208u64, n_hash(17u8));
/// assert_eq!(824633720832u64, n_hash(18u8));
/// assert_eq!(3298534883328u64, n_hash(19u8));
/// assert_eq!(13194139533312u64, n_hash(20u8));
/// assert_eq!(52776558133248u64, n_hash(21u8));
/// assert_eq!(211106232532992u64, n_hash(22u8));
/// assert_eq!(844424930131968u64, n_hash(23u8));
/// assert_eq!(3377699720527872u64, n_hash(24u8));
/// assert_eq!(13510798882111488u64, n_hash(25u8));
/// assert_eq!(54043195528445952u64, n_hash(26u8));
/// assert_eq!(216172782113783808u64, n_hash(27u8));
/// assert_eq!(864691128455135232u64, n_hash(28u8));
/// assert_eq!(3458764513820540928u64, n_hash(29u8));
/// ```
/// 
#[inline]
pub fn n_hash(depth: u8) -> u64 {
  check_depth(depth);
  n_hash_unsafe(depth)
}

/// Same as [n_hash](fn.n_hash.html) except that this version does not panic if the given `depth` is
/// out of range.
#[inline]
pub fn n_hash_unsafe(depth: u8) -> u64 { 12u64 << (depth << 1u8) }

/// Returns `true` if the function [best_starting_depth](fn.best_starting_depth.html) is valid
/// for the given argument `d_max_rad`. So if `d_max_rad < ~48 deg`. `d_max_rad` is given in radians.
/// 
/// ```rust
/// use cdshealpix::{has_best_starting_depth};
/// use std::f64::consts::PI;
/// 
/// assert!(!has_best_starting_depth(PI / 3f64));
/// assert!(has_best_starting_depth(PI / 4f64));
/// ```
#[inline]
pub fn has_best_starting_depth(d_max_rad: f64) -> bool {
    d_max_rad < SMALLER_EDGE2OPEDGE_DIST[0]
}

/// Returns the the smallest depth (in `[0, 29]`) at which a shape having the given largest distance
/// from its center to a border overlaps a maximum of 9 cells (the cell containing the center of
/// the shape plus the 8 neighbouring cells).  
/// Info: internally, unrolled binary search loop on 30 pre-computed values (one by depth).

/// @return -1 if the given distance is very large (> ~48deg), else returns the smallest depth
/// (in [0, 29]) at which a shape having the given largest distance from its center to a border
/// overlaps a maximum of 9 cells (the cell containing the center of the shape plus the 8
/// neighbouring cells).
/// 
/// # Input
/// - `d_max_rad` largest possible distance, in radians, between the center and the border of a shape
///
/// # Output
/// - `depth` = the smallest depth (in `[0, 29]`) at which a shape having the given largest distance 
/// from its center to a border overlaps a maximum of 9 cells (the cell containing the center of the
/// shape plus the 8 neighbouring cells).
///
/// # Panics
/// If the given distance is very large (> ~48deg), this function is not valid since the 12 base
/// cells could be overlaped by the shape 
/// (see [has_best_starting_depth](fn.has_best_starting_depth.html)). Thus it panics.
///
/// # Examples
///
/// ```rust
/// use cdshealpix::{best_starting_depth};
/// use std::f64::consts::PI;
///
/// assert_eq!(0, best_starting_depth(PI / 4f64)); // 45 deg
/// assert_eq!(5, best_starting_depth(0.0174533)); //  1 deg
/// assert_eq!(7, best_starting_depth(0.0043632)); // 15 arcmin
/// assert_eq!(9, best_starting_depth(0.0013));    // 4.469 arcmin
/// assert_eq!(15, best_starting_depth(1.454E-5)); // 3 arcsec
/// assert_eq!(20, best_starting_depth(6.5E-7));   // 0.134 arcsec
/// assert_eq!(22, best_starting_depth(9.537E-8)); // 20 mas
/// ```
#[inline]
pub fn best_starting_depth(d_max_rad: f64) -> u8 { // Could have used an Option
    assert!(d_max_rad < SMALLER_EDGE2OPEDGE_DIST[0],
            "Too large value, use first function has_best_starting_depth");
    // Unrolled binary search loop
    if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[29] {
        29
    } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[15] {
        if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[22] {
            if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[25] {
                if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[27] {
                    if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[28] {
                        28
                    } else {
                        27
                    } 
                } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[26] {
                    26
                } else {
                    25
                } 
            } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[24] {
                24
            } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[23] {
                23
            } else {
                22
            }
        } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[18] {
            if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[20] {
                if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[21] {
                    21
                } else {
                    20
                }
            } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[19] {
                19
            } else {
                18
            }
        } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[17] {
            17
        } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[16] {
            16
        } else {
            15
        }
    } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[7] {
        if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[11] {
            if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[13] {
                if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[14] {
                    14
                } else {
                    13
                }
            } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[12] {
                12
            } else {
                11
            }
        } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[9] {
            if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[10] {
                10
            } else {
                9
            }
        } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[8] {
            8
        } else {
            7
        }
    } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[3] {
        if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[5] {
            if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[6] {
                6
            } else {
                5
            }
        } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[4] {
            4
        } else {
            3
        }
    } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[2] {
        2
    } else if d_max_rad < SMALLER_EDGE2OPEDGE_DIST[1] {
        1
    } else {
        0
    }
}

/// Performs the HEALPix projection: `(x, y) = proj(lon, lat)`.  
/// The chosen scale is such that: base cell vertices and center coordinates are integers;
/// the distance from a cell center to its vertices equals one.  
/// This projection is multi-purpose in the sense that if `lon` is in `[-pi, pi]`, then
/// `x` is in `[-4, 4]` and if `lon` is in `[0, 2pi]`, then `x` is in `[0, 8]`.  
/// It means that a same position on the sphere can lead to different positions in the projected 
/// Euclidean plane.
/// 
/// Simplified projection formulae are:
///  - Equatorial region
/// ```math
/// \boxed{
///   \left\{
///     \begin{array}{lcl}
///       X & = & \alpha \times \frac{4}{\pi} \\
///       Y & = & \sin(\delta) \times \frac{3}{2}
///     \end{array}
///   \right.
/// }
/// \Rightarrow
/// \left\{
///   \begin{array}{lcl}
///     \alpha \in [0, 2\pi] & \leadsto &  X \in [0, 8] \\
///     \sin\delta \in [-\frac{2}{3}, \frac{2}{3}] & \leadsto & Y \in [-1, 1]
///   \end{array}
/// \right.
/// ```
///  - Polar caps:
/// ```math
/// \boxed{
///   \left\{
///     \begin{array}{lcl}
///       t & = & \sqrt{3(1-\sin\delta)} \\
///       X & = & (\alpha\frac{4}{\pi} - 1)t+1 \\
///       Y & = & 2 - t
///     \end{array}
///   \right.
/// }
/// \Rightarrow
/// \left\{
///   \begin{array}{l}
///     \alpha \in [0, \frac{\pi}{2}] \\
///     \sin\delta \in ]\frac{2}{3}, 1]
///   \end{array}
/// \right.
/// \leadsto
/// \begin{array}{l}
///    t \in [0, 1[  \\
///    X \in ]0, 2[ \\
///    Y \in ]1, 2]
/// \end{array}
/// ```
/// 
/// It is the responsibility of the caller to homogenize the result according to its needs.
/// ![Proj](https://raw.githubusercontent.com/cds-astro/cds-healpix-rust/master/resources/4doc/hpx_proj.png)
///
/// # Inputs
/// - `lon` longitude in radians, support positive and negative reasonably large values
///   (naive approach, no Cody-Waite nor Payne Hanek range reduction).
/// - `lat` latitude in radians, must be in `[-pi/2, pi/2]`
///
/// # Output
/// - `(x, y)` the projected planar Euclidean coordinates of the point
///    of given coordinates `(lon, lat)` on the unit sphere
///     - `lon` &le; `0` => `x in [-8, 0]`
///     - `lon` &ge; `0` => `x in [0, 8]`
///     - `y in [-2, 2]`
///
/// # Panics
/// If `lat` **not in** `[-pi/2, pi/2]`, this method panics.
///
/// # Examples
/// To obtain the WCS projection (see Calabretta2007), you can write:
/// ```rust
/// use cdshealpix::{HALF_PI, PI_OVER_FOUR, proj};
/// use std::f64::consts::{PI};
///
/// let lon = 25.1f64;
/// let lat = 46.7f64;
///
/// let (mut x, mut y) = proj(lon.to_radians(), lat.to_radians());
/// if x > 4f64 {
///   x -= 8f64;
/// }
/// x *= PI_OVER_FOUR;
/// y *= PI_OVER_FOUR;
///
/// assert!(-PI <= x && x <= PI);
/// assert!(-HALF_PI <= y && y <= HALF_PI);
/// ```
///
/// Other test example:
/// ```rust
/// use cdshealpix::{TRANSITION_LATITUDE, HALF_PI, PI_OVER_FOUR, proj};
///
/// let (x, y) = proj(0.0, 0.0);
/// assert_eq!(0f64, x);
/// assert_eq!(0f64, y);
///
/// assert_eq!((0.0, 1.0), proj(0.0 * HALF_PI, TRANSITION_LATITUDE));
///
/// fn dist(p1: (f64, f64), p2: (f64, f64)) -> f64 {
///     f64::sqrt((p2.0 - p1.0) * (p2.0 - p1.0) + (p2.1 - p1.1) * (p2.1 - p1.1))
/// }
/// assert!(dist((0.0, 0.0), proj(0.0 * HALF_PI, 0.0)) < 1e-15);
/// assert!(dist((0.0, 1.0), proj(0.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((1.0, 2.0), proj(0.0 * HALF_PI + PI_OVER_FOUR, HALF_PI)) < 1e-15);
/// assert!(dist((2.0, 1.0), proj(1.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((3.0, 2.0), proj(1.0 * HALF_PI + PI_OVER_FOUR, HALF_PI)) < 1e-15);
/// assert!(dist((4.0, 1.0), proj(2.0 * HALF_PI , TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((5.0, 2.0), proj(2.0 * HALF_PI + PI_OVER_FOUR, HALF_PI)) < 1e-15);
/// assert!(dist((6.0, 1.0), proj(3.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((7.0, 2.0), proj(3.0 * HALF_PI + PI_OVER_FOUR, HALF_PI)) < 1e-15);
/// assert!(dist((0.0, 1.0), proj(4.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((0.0, 0.0), proj(4.0 * HALF_PI, 0.0)) < 1e-15);
/// assert!(dist((0.0, -1.0), proj(0.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((1.0, -2.0), proj(0.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI)) < 1e-15);
/// assert!(dist((2.0, -1.0), proj(1.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((3.0, -2.0), proj(1.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI)) < 1e-15);
/// assert!(dist((4.0, -1.0), proj(2.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((5.0, -2.0), proj(2.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI)) < 1e-15);
/// assert!(dist((6.0, -1.0), proj(3.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((7.0, -2.0), proj(3.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI)) < 1e-15);
/// assert!(dist((0.0, -1.0), proj(4.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
///
/// assert!(dist((-0.0, 0.0), proj(-0.0 * HALF_PI, 0.0)) < 1e-15);
/// assert!(dist((-0.0, 1.0), proj(-0.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((-1.0, 2.0), proj(-0.0 * HALF_PI - PI_OVER_FOUR, HALF_PI)) < 1e-15);
/// assert!(dist((-2.0, 1.0), proj(-1.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((-3.0, 2.0), proj(-1.0 * HALF_PI - PI_OVER_FOUR, HALF_PI)) < 1e-15);
/// assert!(dist((-4.0, 1.0), proj(-2.0 * HALF_PI , TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((-5.0, 2.0), proj(-2.0 * HALF_PI - PI_OVER_FOUR, HALF_PI)) < 1e-15);
/// assert!(dist((-6.0, 1.0), proj(-3.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((-7.0, 2.0), proj(-3.0 * HALF_PI - PI_OVER_FOUR, HALF_PI)) < 1e-15);
/// assert!(dist((-0.0, 1.0), proj(-4.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((-0.0, 0.0), proj(-4.0 * HALF_PI, 0.0)) < 1e-15);
/// assert!(dist((-0.0, -1.0), proj(-0.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((-1.0, -2.0), proj(-0.0 * HALF_PI - PI_OVER_FOUR, -HALF_PI)) < 1e-15);
/// assert!(dist((-2.0, -1.0), proj(-1.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((-3.0, -2.0), proj(-1.0 * HALF_PI - PI_OVER_FOUR, -HALF_PI)) < 1e-15);
/// assert!(dist((-4.0, -1.0), proj(-2.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((-5.0, -2.0), proj(-2.0 * HALF_PI - PI_OVER_FOUR, -HALF_PI)) < 1e-15);
/// assert!(dist((-6.0, -1.0), proj(-3.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
/// assert!(dist((-7.0, -2.0), proj(-3.0 * HALF_PI - PI_OVER_FOUR, -HALF_PI)) < 1e-15);
/// assert!(dist((-0.0, -1.0), proj(-4.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
/// ```
#[inline]
pub fn proj(lon: f64, lat: f64) -> (f64, f64) {
    check_lat(lat);
    let lon = abs_sign_decompose(lon);
    let lat = abs_sign_decompose(lat);
    let x = pm1_offset_decompose(lon.abs * FOUR_OVER_PI);
    let mut xy = (x.pm1, lat.abs);
    if is_in_equatorial_region(lat.abs) {
        proj_cea(&mut xy);
    } else {
        proj_collignon(&mut xy);
    }
    apply_offset_and_signs(&mut xy, x.offset, lon.sign, lat.sign);
    xy
}


/// Unproject the given HEALPix projected points.  
/// This unprojection is multi-purpose in the sense that:
///  - if input `x` in `[-8, 0[`, then output `lon` in `[-2pi, 0]`
///  - if input `x` in `[ 0, 8]`, then output `lon` in `[0, 2pi]`
///  - output `lat` always in `[-pi/2, pi/2]`
/// 
/// # Inputs
///  - `x` the projected coordinate along the x-axis, supports positive and negative reasonably 
///        large values with a naive approach (no Cody-Waite nor Payne Hanek range reduction).
///  - `y` the projected coordinate along te x-axis, must be in `[-2, 2]`
/// 
/// # Output
/// -  `(lon, lat)` in radians, the position on the unit sphere whose projected coordinates are 
///    the input coordinates `(x, y)`.  
///   - if `x <= 0`, then `lon` in `[-2pi, 0]`;
///   - else if `x >= 0`, the  `lon` in `[0, 2pi]`
///   - `lat` always in `[-pi/2, pi/2]`.
///
/// # Panics
/// If `y` **not in** `[-2, 2]`, this method panics.
///
/// # Examples
/// To obtain the WCS un-projection (see Calabretta2007), you can write:
/// ```rust
/// use cdshealpix::{HALF_PI, FOUR_OVER_PI, unproj};
/// use std::f64::consts::{PI};
///
/// let x = 2.1f64;
/// let y = 0.36f64;
///
/// let (mut lon, mut lat) = unproj(x * FOUR_OVER_PI, y * FOUR_OVER_PI);
/// if lon < 0f64 {
///     lon += 2f64 * PI;
/// }
///
/// assert!(0f64 <= lon && lon <= 2f64 * PI);
/// assert!(-HALF_PI <= lat && lat <= HALF_PI);
/// ```
///
/// Other test example:
/// ```rust
/// use cdshealpix::{TRANSITION_LATITUDE, HALF_PI, PI_OVER_FOUR, proj, unproj};
///
/// fn dist(p1: (f64, f64), p2: (f64, f64)) -> f64 {
///     let sindlon = f64::sin(0.5 * (p2.0 - p1.0));
///     let sindlat = f64::sin(0.5 * (p2.1 - p1.1));
///     2f64 * f64::asin(f64::sqrt(sindlat * sindlat + p1.1.cos() * p2.1.cos() * sindlon * sindlon))
/// }
/// 
/// let points: [(f64, f64); 40] = [
///     (0.0 * HALF_PI, 0.0),
///     (1.0 * HALF_PI, TRANSITION_LATITUDE),
///     (2.0 * HALF_PI, TRANSITION_LATITUDE),
///     (3.0 * HALF_PI, TRANSITION_LATITUDE),
///     (4.0 * HALF_PI, TRANSITION_LATITUDE),
///     (0.0 * HALF_PI + PI_OVER_FOUR, HALF_PI),
///     (1.0 * HALF_PI + PI_OVER_FOUR, HALF_PI),
///     (2.0 * HALF_PI + PI_OVER_FOUR, HALF_PI),
///     (3.0 * HALF_PI + PI_OVER_FOUR, HALF_PI),
///     (4.0 * HALF_PI + PI_OVER_FOUR, HALF_PI),
///     (0.0 * HALF_PI, 0.0),
///     (1.0 * HALF_PI, -TRANSITION_LATITUDE),
///     (2.0 * HALF_PI, -TRANSITION_LATITUDE),
///     (3.0 * HALF_PI, -TRANSITION_LATITUDE),
///     (4.0 * HALF_PI, -TRANSITION_LATITUDE),
///     (0.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI),
///     (1.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI),
///     (2.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI),
///     (3.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI),
///     (4.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI),
///     (-0.0 * HALF_PI, 0.0),
///     (-1.0 * HALF_PI, TRANSITION_LATITUDE),
///     (-2.0 * HALF_PI, TRANSITION_LATITUDE),
///     (-3.0 * HALF_PI, TRANSITION_LATITUDE),
///     (-4.0 * HALF_PI, TRANSITION_LATITUDE),
///     (-0.0 * HALF_PI + PI_OVER_FOUR, HALF_PI),
///     (-1.0 * HALF_PI + PI_OVER_FOUR, HALF_PI),
///     (-2.0 * HALF_PI + PI_OVER_FOUR, HALF_PI),
///     (-3.0 * HALF_PI + PI_OVER_FOUR, HALF_PI),
///     (-4.0 * HALF_PI + PI_OVER_FOUR, HALF_PI),
///     (-0.0 * HALF_PI, 0.0),
///     (-1.0 * HALF_PI, -TRANSITION_LATITUDE),
///     (-2.0 * HALF_PI, -TRANSITION_LATITUDE),
///     (-3.0 * HALF_PI, -TRANSITION_LATITUDE),
///     (-4.0 * HALF_PI, -TRANSITION_LATITUDE),
///     (-0.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI),
///     (-1.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI),
///     (-2.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI),
///     (-3.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI),
///     (-4.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI)
/// ];
/// 
/// for (lon, lat) in points.iter() {
///     let (x, y): (f64, f64) = proj(*lon, *lat);
///     assert!(dist((*lon, *lat), unproj(x, y)) < 1e-15);
/// }
/// ```
#[inline]
pub fn unproj(x: f64, y: f64) -> (f64, f64) {
    check_y(y);
    let x = abs_sign_decompose(x);
    let y = abs_sign_decompose(y);
    let lon = pm1_offset_decompose(x.abs);
    let mut lonlat= (lon.pm1, y.abs);
    if is_in_projected_equatorial_region(y.abs) {
        deproj_cea(&mut lonlat);
    } else {
        deproj_collignon(&mut lonlat);
    }
    apply_offset_and_signs(&mut lonlat, lon.offset, x.sign, y.sign);
    lonlat.0 *= PI_OVER_FOUR;
    lonlat
}


/// Verify that the latitude is in [-PI/2, PI/2], panics if not.
#[inline]
fn check_lat(lat: f64) {
    assert!(-HALF_PI <= lat && lat <= HALF_PI);
}

/// Verify that the projected y coordinate is in [-2, 2], panics if not.
#[inline]
fn check_y(y: f64) { assert!(-2f64 <= y && y <= 2f64); }

/// Returns `true` if the point of given (absolute value of) latitude is in the equatorial region,
/// and `false` if it is located in one of the two polar caps
#[inline]
pub fn is_in_equatorial_region(abs_lat: f64) -> bool {
    abs_lat <= TRANSITION_LATITUDE
}

/// Returns `true` if the point of given (absolute value of) y coordinate in the projected plane
/// is in the equatorial region, and `false` if it is located in one of the two polar caps
#[inline]
pub fn is_in_projected_equatorial_region(abs_y: f64) -> bool { abs_y <= 1f64 }

// Returns the absolute value of the given double together with its bit of sign
struct AbsAndSign {
    abs: f64,
    sign: u64,
}
#[inline]
fn abs_sign_decompose(x: f64) -> AbsAndSign {
    let bits = f64::to_bits(x);
    AbsAndSign {
        abs: f64::from_bits(bits & F64_BUT_SIGN_BIT_MASK),
        sign: bits & F64_SIGN_BIT_MASK,
    }
}

// Decompose the given positive real value in
// --* an integer offset in [1, 3, 5, 7] (*PI/4) and
// --* a real value in [-1.0, 1.0] (*PI/4)
struct OffsetAndPM1 {
    offset: u8, // = 1, 3, 5 or 7
    pm1: f64,   // in [-1.0, 1.0]
}
#[inline]
fn pm1_offset_decompose(x: f64) -> OffsetAndPM1 {
    let floor: u8 = x as u8;
    let odd_floor: u8 = floor | 1u8;
    OffsetAndPM1 {
        offset: odd_floor & 7u8, // value modulo 7
        pm1: x - (odd_floor as f64),
    }
}

// Cylindrical Equal Area projection
#[inline]
fn proj_cea(xy: &mut (f64, f64)) {
    let (_, ref mut y) = *xy;
    *y = f64::sin(*y) * ONE_OVER_TRANSITION_Z;
}
#[inline]
fn deproj_cea(lonlat: &mut (f64, f64)) {
    let (_, ref mut lat) = *lonlat;
    *lat = f64::asin((*lat) * TRANSITION_Z);
}

// Collignon projection
#[inline]
fn proj_collignon(xy: &mut (f64, f64)) {
    let (ref mut x, ref mut y) = *xy;
    *y = SQRT6 * f64::cos(HALF * *y + PI_OVER_FOUR);
    *x *= *y;
    *y = 2f64 - *y;
}
#[inline]
fn deproj_collignon(lonlat: &mut (f64, f64)) {
    let (ref mut lon, ref mut lat) = *lonlat;
    *lat = 2f64 - *lat;
    if is_not_near_from_pole(*lat) { // Rare, so few risks of branch miss-prediction
        *lon /= *lat;
        deal_with_numerical_approx_in_edges(lon);
    } // in case of pole, lon = lat = 0 (we avoid NaN due to division by lat=0)
    *lat *= ONE_OVER_SQRT6;
    *lat = 2f64 * f64::acos(*lat) - HALF_PI;
}

#[inline]
fn is_not_near_from_pole(sqrt_of_three_time_one_minus_sin_of: f64) -> bool {
    // In case of pole: x = y = 0
    sqrt_of_three_time_one_minus_sin_of > EPS_POLE
}

#[inline]
fn deal_with_numerical_approx_in_edges(lon: &mut f64) {
    if *lon > 1f64 {
        *lon = 1f64;
    } else if *lon < -1f64 {
        *lon = -1f64;
    }
}

// Shift x by the given offset and apply lon and lat signs to x and y respectively
#[inline]
fn apply_offset_and_signs(ab: &mut (f64, f64), off: u8, a_sign: u64, b_sign: u64) {
    let (ref mut a, ref mut b) = *ab;
    *a += off as f64;
    *a = f64::from_bits(f64::to_bits(*a) | a_sign);
    *b = f64::from_bits(f64::to_bits(*b) | b_sign);
}

// Import module compass point
pub mod compass_point;
use crate::compass_point::{MainWind};
use crate::compass_point::MainWind::*;

/// Compute the base cell value which is the neighbour of the given base cell, in the given direction.  
/// There is no neighbour:
/// - in the North and South directions for the equatorial region cells (i.e. cells 4, 5, 6 and 7)
/// - in the East and West directions for:
///   - the north polar cap cells (i.e. cells 0, 1, 2 and 3)
///   - the south polar cap cells (i.e. cells 8, 9, 10 and 11)  
pub fn neighbour(base_cell: u8, direction: MainWind) -> Option<u8> {
  if direction == MainWind::C {
    Some(base_cell)
  } else {
    let d0h_mod_4 = base_cell & 3_u8;  // <=> base_cell modulo 4
    match base_cell >> 2 { // <=> basce_cell / 4
      0 => ncp_neighbour(d0h_mod_4, direction),
      1 => eqr_neighbour(d0h_mod_4, direction),
      2 => spc_neighbour(d0h_mod_4, direction),
      _ => panic!("Base cell must be in [0, 12["),
    }
  }
}

fn ncp_neighbour(d0h_mod_4: u8, direction: MainWind) -> Option<u8> {
  match direction {
     S => base_cell_opt(iden(d0h_mod_4), 2),
    SE => base_cell_opt(next(d0h_mod_4), 1),
    SW => base_cell_opt(iden(d0h_mod_4), 1),
    NE => base_cell_opt(next(d0h_mod_4), 0),
    NW => base_cell_opt(prev(d0h_mod_4), 0),
     N => base_cell_opt(oppo(d0h_mod_4), 0),
    _ => None,
  }  
}

fn eqr_neighbour(d0h_mod_4: u8, direction: MainWind) -> Option<u8> {
  match direction {
    SE => base_cell_opt(iden(d0h_mod_4), 2),
     E => base_cell_opt(next(d0h_mod_4), 1),
    SW => base_cell_opt(prev(d0h_mod_4), 2),
    NE => base_cell_opt(iden(d0h_mod_4), 0),
     W => base_cell_opt(prev(d0h_mod_4), 1),
    NW => base_cell_opt(prev(d0h_mod_4), 0),
    _ => None,
  }
}

fn spc_neighbour(d0h_mod_4: u8, direction: MainWind) -> Option<u8> {
  match direction {
     S => base_cell_opt(oppo(d0h_mod_4), 2),
    SE => base_cell_opt(next(d0h_mod_4), 2),
    SW => base_cell_opt(prev(d0h_mod_4), 2),
    NE => base_cell_opt(next(d0h_mod_4), 1),
    NW => base_cell_opt(iden(d0h_mod_4), 1),
     N => base_cell_opt(iden(d0h_mod_4), 0),
    _ => None,
  }
}

/// Returns (mod4 - 1) in [0, 2], and 3 if mod4 == 0 (i.e. the previous value in [0, 3] range)
#[inline]
fn prev(mod4: u8) -> u8 {
  debug_assert!(mod4 < 4);
  (((mod4 as i8) - 1) & 3) as u8
}

/// Returns (mod4 + 1) in [1, 3], and 0 if mod4 == 3 (i.e. the next value in [0, 3] range)
#[inline]
fn next(mod4: u8) -> u8 {
  debug_assert!(mod4 < 4);
  (mod4 + 1) & 3
}

/// Returns (mod4 + 2) in [2, 3], and 0 if mod4 == 2 and 1 if mod4 == 3 (i.e. the opposite value in [0, 3] range)
#[inline]
fn oppo(mod4: u8) -> u8 {
  debug_assert!(mod4 < 4);
  (mod4 + 2) & 3
}

/// Returns the input value: useless, just used to improve code legibility
#[inline]
fn iden(mod4: u8) -> u8 {
  debug_assert!(mod4 < 4);
  mod4
}

#[inline]
fn base_cell_opt(i: u8, j: u8) -> Option<u8> {
  Some(base_cell(i, j))
}

/// Compute the base cell from its (i, j) coordinates:
/// - i: index along the longitude axis ( = base_cell modulo 4)
/// - j: index along the latitude axis ( = base_cell / 4)
///   - = 0 for the cells covering the north polar cap
///   - = 1 for the cells with are only in the equatorial region
///   - = 2 for the cells covering the south polar cap
#[inline]
fn base_cell(i: u8, j: u8) -> u8 {
  debug_assert!(i < 4 && j < 3);
  (j << 2) + i
}


/// Module containing NESTED scheme methods
pub mod nested;

/// No need to make those public!
mod sph_geom;
mod special_points_finder;

#[cfg(test)]
mod tests {
  use super::*;

  #[test]
  fn testok_nside() {
    assert_eq!(1, nside(0));
    assert_eq!(2, nside(1));
    assert_eq!(4, nside(2));
    assert_eq!(8, nside(3));
    assert_eq!(16, nside(4));
    assert_eq!(32, nside(5));
    assert_eq!(64, nside(6));
    assert_eq!(128, nside(7));
    assert_eq!(256, nside(8));
    assert_eq!(512, nside(9));
    assert_eq!(1024, nside(10));
    assert_eq!(2048, nside(11));
    assert_eq!(4096, nside(12));
    assert_eq!(8192, nside(13));
    assert_eq!(16384, nside(14));
    assert_eq!(32768, nside(15));
    assert_eq!(65536, nside(16));
    assert_eq!(131072, nside(17));
    assert_eq!(262144, nside(18));
    assert_eq!(524288, nside(19));
    assert_eq!(1048576, nside(20));
    assert_eq!(2097152, nside(21));
    assert_eq!(4194304, nside(22));
    assert_eq!(8388608, nside(23));
    assert_eq!(16777216, nside(24));
    assert_eq!(33554432, nside(25));
    assert_eq!(67108864, nside(26));
    assert_eq!(134217728, nside(27));
    assert_eq!(268435456, nside(28));
    assert_eq!(536870912, nside(29));
  }

  #[test]
  #[should_panic]
  fn testpanic_nside() {
    nside(30);
  }

  #[test]
  fn testok_proj() {
    fn dist(p1: (f64, f64), p2: (f64, f64)) -> f64 {
      f64::sqrt((p2.0 - p1.0) * (p2.0 - p1.0) + (p2.1 - p1.1) * (p2.1 - p1.1))
    }
    // println!("{:?}", proj(0.0, TRANSITION_LATITUDE));
    // println!("{}", dist((0.0, 1.0), proj(0.0, TRANSITION_LATITUDE)));
    assert!(dist((0.0, 0.0), proj(0.0 * HALF_PI, 0.0)) < 1e-15);
    assert!(dist((0.0, 1.0), proj(0.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((1.0, 2.0), proj(0.0 * HALF_PI + PI_OVER_FOUR, HALF_PI)) < 1e-15);
    assert!(dist((2.0, 1.0), proj(1.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((3.0, 2.0), proj(1.0 * HALF_PI + PI_OVER_FOUR, HALF_PI)) < 1e-15);
    assert!(dist((4.0, 1.0), proj(2.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((5.0, 2.0), proj(2.0 * HALF_PI + PI_OVER_FOUR, HALF_PI)) < 1e-15);
    assert!(dist((6.0, 1.0), proj(3.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((7.0, 2.0), proj(3.0 * HALF_PI + PI_OVER_FOUR, HALF_PI)) < 1e-15);
    assert!(dist((0.0, 1.0), proj(4.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((0.0, 0.0), proj(4.0 * HALF_PI, 0.0)) < 1e-15);
    assert!(dist((0.0, -1.0), proj(0.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((1.0, -2.0), proj(0.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI)) < 1e-15);
    assert!(dist((2.0, -1.0), proj(1.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((3.0, -2.0), proj(1.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI)) < 1e-15);
    assert!(dist((4.0, -1.0), proj(2.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((5.0, -2.0), proj(2.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI)) < 1e-15);
    assert!(dist((6.0, -1.0), proj(3.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((7.0, -2.0), proj(3.0 * HALF_PI + PI_OVER_FOUR, -HALF_PI)) < 1e-15);
    assert!(dist((0.0, -1.0), proj(4.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);

    assert!(dist((-0.0, 0.0), proj(-0.0 * HALF_PI, 0.0)) < 1e-15);
    assert!(dist((-0.0, 1.0), proj(-0.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((-1.0, 2.0), proj(-0.0 * HALF_PI - PI_OVER_FOUR, HALF_PI)) < 1e-15);
    assert!(dist((-2.0, 1.0), proj(-1.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((-3.0, 2.0), proj(-1.0 * HALF_PI - PI_OVER_FOUR, HALF_PI)) < 1e-15);
    assert!(dist((-4.0, 1.0), proj(-2.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((-5.0, 2.0), proj(-2.0 * HALF_PI - PI_OVER_FOUR, HALF_PI)) < 1e-15);
    assert!(dist((-6.0, 1.0), proj(-3.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((-7.0, 2.0), proj(-3.0 * HALF_PI - PI_OVER_FOUR, HALF_PI)) < 1e-15);
    assert!(dist((-0.0, 1.0), proj(-4.0 * HALF_PI, TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((-0.0, 0.0), proj(-4.0 * HALF_PI, 0.0)) < 1e-15);
    assert!(dist((-0.0, -1.0), proj(-0.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((-1.0, -2.0), proj(-0.0 * HALF_PI - PI_OVER_FOUR, -HALF_PI)) < 1e-15);
    assert!(dist((-2.0, -1.0), proj(-1.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((-3.0, -2.0), proj(-1.0 * HALF_PI - PI_OVER_FOUR, -HALF_PI)) < 1e-15);
    assert!(dist((-4.0, -1.0), proj(-2.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((-5.0, -2.0), proj(-2.0 * HALF_PI - PI_OVER_FOUR, -HALF_PI)) < 1e-15);
    assert!(dist((-6.0, -1.0), proj(-3.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
    assert!(dist((-7.0, -2.0), proj(-3.0 * HALF_PI - PI_OVER_FOUR, -HALF_PI)) < 1e-15);
    assert!(dist((-0.0, -1.0), proj(-4.0 * HALF_PI, -TRANSITION_LATITUDE)) < 1e-15);
  }

  #[test]
  #[should_panic]
  fn testpanic_proj_1() {
    proj(0.0, -1.58);
  }

  #[test]
  #[should_panic]
  fn testpanic_proj_2() {
    proj(3.14159, -1.58);
  }

  #[test]
  fn testok_shs() {
    let ang_dist = 0.24;
    let shs = to_squared_half_segment(ang_dist);
    let ang_dist_2 = sphe_dist(shs);
    assert!((ang_dist - ang_dist_2).abs() < 1e-4);
  }
  
}