Source

yt / yt / visualization / volume_rendering / camera.py

Full commit
   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
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
"""
Import the components of the volume rendering extension

Author: Matthew Turk <matthewturk@gmail.com>
Affiliation: KIPAC/SLAC/Stanford
Homepage: http://yt-project.org/
License:
  Copyright (C) 2009 Matthew Turk.  All Rights Reserved.

  This file is part of yt.

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

  This program is distributed in the hope that it will be useful,
  but WITHOUT ANY WARRANTY; without even the implied warranty of
  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
  GNU General Public License for more details.

  You should have received a copy of the GNU General Public License
  along with this program.  If not, see <http://www.gnu.org/licenses/>.
"""

import __builtin__
import numpy as np

from yt.funcs import *
from yt.utilities.math_utils import *

from .grid_partitioner import HomogenizedVolume
from .transfer_functions import ProjectionTransferFunction

from yt.utilities.lib import \
    arr_vec2pix_nest, arr_pix2vec_nest, \
    arr_ang2pix_nest, arr_fisheye_vectors
from yt.utilities.math_utils import get_rotation_matrix
from yt.utilities.orientation import Orientation
from yt.data_objects.api import ImageArray
from yt.visualization.image_writer import write_bitmap, write_image
from yt.data_objects.data_containers import data_object_registry
from yt.utilities.parallel_tools.parallel_analysis_interface import \
    ParallelAnalysisInterface, ProcessorPool
from yt.utilities.amr_kdtree.api import AMRKDTree

from yt.utilities.lib import \
    PartitionedGrid, ProjectionSampler, VolumeRenderSampler, \
    LightSourceRenderSampler, InterpolatedProjectionSampler, \
    arr_vec2pix_nest, arr_pix2vec_nest, arr_ang2pix_nest, \
    pixelize_healpix, arr_fisheye_vectors

class Camera(ParallelAnalysisInterface):

    _sampler_object = VolumeRenderSampler

    def __init__(self, center, normal_vector, width,
                 resolution, transfer_function,
                 north_vector = None, steady_north=False,
                 volume = None, fields = None,
                 log_fields = None,
                 sub_samples = 5, pf = None,
                 use_kd=True, l_max=None, no_ghost=True,
                 tree_type='domain',
                 le=None, re=None, use_light=False):
        r"""A viewpoint into a volume, for volume rendering.

        The camera represents the eye of an observer, which will be used to
        generate ray-cast volume renderings of the domain.

        Parameters
        ----------
        center : array_like
            The current "center" of the view port -- the focal point for the
            camera.
        normal_vector : array_like
            The vector between the camera position and the center.
        width : float or list of floats
            The current width of the image.  If a single float, the volume is
            cubical, but if not, it is left/right, top/bottom, front/back.
        resolution : int or list of ints
            The number of pixels in each direction.
        north_vector : array_like, optional
            The 'up' direction for the plane of rays.  If not specific, calculated
            automatically.
        steady_north : bool, optional
            Boolean to control whether to normalize the north_vector
            by subtracting off the dot product of it and the normal
            vector.  Makes it easier to do rotations along a single
            axis.  If north_vector is specified, is switched to
            True. Default: False
        volume : `yt.extensions.volume_rendering.HomogenizedVolume`, optional
            The volume to ray cast through.  Can be specified for finer-grained
            control, but otherwise will be automatically generated.
        fields : list of fields, optional
            This is the list of fields we want to volume render; defaults to
            Density.
        log_fields : list of bool, optional
            Whether we should take the log of the fields before supplying them to
            the volume rendering mechanism.
        sub_samples : int, optional
            The number of samples to take inside every cell per ray.
        pf : `~yt.data_objects.api.StaticOutput`
            For now, this is a require parameter!  But in the future it will become
            optional.  This is the parameter file to volume render.
        use_kd: bool, optional
            Specifies whether or not to use a kd-Tree framework for
            the Homogenized Volume and ray-casting.  Default to True.
        l_max: int, optional
            Specifies the maximum level to be rendered.  Also
            specifies the maximum level used in the kd-Tree
            construction.  Defaults to None (all levels), and only
            applies if use_kd=True.
        no_ghost: bool, optional
            Optimization option.  If True, homogenized bricks will
            extrapolate out from grid instead of interpolating from
            ghost zones that have to first be calculated.  This can
            lead to large speed improvements, but at a loss of
            accuracy/smoothness in resulting image.  The effects are
            less notable when the transfer function is smooth and
            broad. Default: False
        tree_type: string, optional
            Specifies the type of kd-Tree to be constructed/cast.
            There are three options, the default being 'domain'. Only
            affects parallel rendering.  'domain' is suggested.

            'domain' - Tree construction/casting is load balanced by
            splitting up the domain into the first N subtrees among N
            processors (N must be a power of 2).  Casting then
            proceeds with each processor rendering their subvolume,
            and final image is composited on the root processor.  The
            kd-Tree is never combined, reducing communication and
            memory overhead. The viewpoint can be changed without
            communication or re-partitioning of the data, making it
            ideal for rotations/spins.

            'breadth' - kd-Tree is first constructed as in 'domain',
            but then combined among all the subtrees.  Rendering is
            then split among N processors (again a power of 2), based
            on the N most expensive branches of the tree.  As in
            'domain', viewpoint can be changed without re-partitioning
            or communication.

            'depth' - kd-Tree is first constructed as in 'domain', but
            then combined among all subtrees.  Rendering is then load
            balanced in a back-to-front manner, splitting up the cost
            as evenly as possible.  If the viewpoint changes,
            additional data might have to be partitioned.  Is also
            prone to longer data IO times.  If all the data can fit in
            memory on each cpu, this can be the fastest option for
            multiple ray casts on the same dataset.
        le: array_like, optional
            Specifies the left edge of the volume to be rendered.
            Currently only works with use_kd=True.
        re: array_like, optional
            Specifies the right edge of the volume to be rendered.
            Currently only works with use_kd=True.

        Examples
        --------

        >>> cam = vr.Camera(c, L, W, (N,N), transfer_function = tf, pf = pf)
        >>> image = cam.snapshot()

        >>> from yt.mods import *
        >>> import yt.visualization.volume_rendering.api as vr
        
        >>> pf = EnzoStaticOutput('DD1701') # Load pf
        >>> c = [0.5]*3 # Center
        >>> L = [1.0,1.0,1.0] # Viewpoint
        >>> W = np.sqrt(3) # Width
        >>> N = 1024 # Pixels (1024^2)

        # Get density min, max
        >>> mi, ma = pf.h.all_data().quantities['Extrema']('Density')[0]
        >>> mi, ma = np.log10(mi), np.log10(ma)

        # Construct transfer function
        >>> tf = vr.ColorTransferFunction((mi-2, ma+2))
        # Sample transfer function with 5 gaussians.  Use new col_bounds keyword.
        >>> tf.add_layers(5,w=0.05, col_bounds = (mi+1,ma), colormap='spectral')
        
        # Create the camera object
        >>> cam = vr.Camera(c, L, W, (N,N), transfer_function=tf, pf=pf) 
        
        # Ray cast, and save the image.
        >>> image = cam.snapshot(fn='my_rendering.png')

        """
        ParallelAnalysisInterface.__init__(self)
        if pf is not None: self.pf = pf
        if not iterable(resolution):
            resolution = (resolution, resolution)
        self.resolution = resolution
        self.sub_samples = sub_samples
        if not iterable(width):
            width = (width, width, width) # left/right, top/bottom, front/back 
        self.orienter = Orientation(normal_vector, north_vector=north_vector, steady_north=steady_north)
        self.rotation_vector = self.orienter.unit_vectors[1]
        self._setup_box_properties(width, center, self.orienter.unit_vectors)
        if fields is None: fields = ["Density"]
        self.fields = fields
        if transfer_function is None:
            transfer_function = ProjectionTransferFunction()
        self.transfer_function = transfer_function
        self.log_fields = log_fields
        self.use_kd = use_kd
        self.l_max = l_max
        self.no_ghost = no_ghost
        self.use_light = use_light
        self.light_dir = None
        self.light_rgba = None
        if self.no_ghost:
            mylog.info('Warning: no_ghost is currently True (default). This may lead to artifacts at grid boundaries.')
        self.tree_type = tree_type
        if volume is None:
            if self.use_kd:
                volume = AMRKDTree(self.pf, l_max=l_max, fields=self.fields, no_ghost=no_ghost, tree_type=tree_type,
                                   log_fields = log_fields, le=le, re=re)
            else:
                volume = HomogenizedVolume(fields, pf = self.pf,
                                           log_fields = log_fields)
        else:
            self.use_kd = isinstance(volume, AMRKDTree)
        self.volume = volume

    def _setup_box_properties(self, width, center, unit_vectors):
        self.width = width
        self.center = center
        self.box_vectors = np.array([unit_vectors[0]*width[0],
                                     unit_vectors[1]*width[1],
                                     unit_vectors[2]*width[2]])
        self.origin = center - 0.5*np.dot(width,unit_vectors)
        self.back_center =  center - 0.5*width[2]*unit_vectors[2]
        self.front_center = center + 0.5*width[2]*unit_vectors[2]         

    def update_view_from_matrix(self, mat):
        pass

    def look_at(self, new_center, north_vector = None):
        r"""Change the view direction based on a new focal point.

        This will recalculate all the necessary vectors and vector planes to orient
        the image plane so that it points at a new location.

        Parameters
        ----------
        new_center : array_like
            The new "center" of the view port -- the focal point for the
            camera.
        north_vector : array_like, optional
            The "up" direction for the plane of rays.  If not specific,
            calculated automatically.
        """
        normal_vector = self.front_center - new_center
        self.orienter.switch_orientation(normal_vector=normal_vector,
                                         north_vector = north_vector)

    def switch_view(self, normal_vector=None, width=None, center=None, north_vector=None):
        r"""Change the view based on any of the view parameters.

        This will recalculate the orientation and width based on any of
        normal_vector, width, center, and north_vector.

        Parameters
        ----------
        normal_vector: array_like, optional
            The new looking vector.
        width: float or array of floats, optional
            The new width.  Can be a single value W -> [W,W,W] or an
            array [W1, W2, W3] (left/right, top/bottom, front/back)
        center: array_like, optional
            Specifies the new center.
        north_vector : array_like, optional
            The 'up' direction for the plane of rays.  If not specific,
            calculated automatically.
        """
        if width is None:
            width = self.width
        if not iterable(width):
            width = (width, width, width) # left/right, tom/bottom, front/back 
        self.width = width
        if center is not None:
            self.center = center
        if north_vector is None:
            north_vector = self.orienter.unit_vectors[1]
        if normal_vector is None:
            normal_vector = self.orienter.normal_vector
        self.orienter.switch_orientation(normal_vector = normal_vector,
                                         north_vector = north_vector)
        self._setup_box_properties(width, self.center, self.orienter.unit_vectors)
    def new_image(self):
        image = np.zeros((self.resolution[0], self.resolution[1], 3), dtype='float64', order='C')
        return image

    def get_sampler_args(self, image):
        rotp = np.concatenate([self.orienter.inv_mat.ravel('F'), self.back_center.ravel()])
        args = (rotp, self.box_vectors[2], self.back_center,
                (-self.width[0]/2.0, self.width[0]/2.0,
                 -self.width[1]/2.0, self.width[1]/2.0),
                image, self.orienter.unit_vectors[0], self.orienter.unit_vectors[1],
                np.array(self.width), self.transfer_function, self.sub_samples)
        return args

    star_trees = None
    def get_sampler(self, args):
        kwargs = {}
        if self.star_trees is not None:
            kwargs = {'star_list': self.star_trees}
        if self.use_light:
            if self.light_dir is None:
                self.set_default_light_dir()
            temp_dir = np.empty(3,dtype='float64')
            temp_dir = self.light_dir[0] * self.orienter.unit_vectors[1] + \
                    self.light_dir[1] * self.orienter.unit_vectors[2] + \
                    self.light_dir[2] * self.orienter.unit_vectors[0]
            if self.light_rgba is None:
                self.set_default_light_rgba()
            sampler = LightSourceRenderSampler(*args, light_dir=temp_dir,
                    light_rgba=self.light_rgba, **kwargs)
        else:
            sampler = self._sampler_object(*args, **kwargs)
        print sampler, kwargs
        return sampler

    def finalize_image(self, image):
        pass

    def _render(self, double_check, num_threads, image, sampler):
        pbar = get_pbar("Ray casting", (self.volume.brick_dimensions + 1).prod(axis=-1).sum())
        total_cells = 0
        if double_check:
            for brick in self.volume.bricks:
                for data in brick.my_data:
                    if np.any(np.isnan(data)):
                        raise RuntimeError

        view_pos = self.front_center + self.orienter.unit_vectors[2] * 1.0e6 * self.width[2]
        for brick in self.volume.traverse(view_pos, self.front_center, image):
            sampler(brick, num_threads=num_threads)
            total_cells += np.prod(brick.my_data[0].shape)
            pbar.update(total_cells)

        pbar.finish()
        image = sampler.aimage
        self.finalize_image(image)
        return image

    def save_image(self, fn, clip_ratio, image):
        if self.comm.rank is 0 and fn is not None:
            image.write_png(fn, clip_ratio=clip_ratio)

    def initialize_source(self):
        return self.volume.initialize_source()

    def get_information(self):
        info_dict = {'fields':self.fields,
                     'type':self.__class__.__name__,
                     'east_vector':self.orienter.unit_vectors[0],
                     'north_vector':self.orienter.unit_vectors[1],
                     'normal_vector':self.orienter.unit_vectors[2],
                     'width':self.width,
                     'dataset':self.pf.fullpath}
        return info_dict

    def snapshot(self, fn = None, clip_ratio = None, double_check = False,
                 num_threads = 0):
        r"""Ray-cast the camera.

        This method instructs the camera to take a snapshot -- i.e., call the ray
        caster -- based on its current settings.

        Parameters
        ----------
        fn : string, optional
            If supplied, the image will be saved out to this before being
            returned.  Scaling will be to the maximum value.
        clip_ratio : float, optional
            If supplied, the 'max_val' argument to write_bitmap will be handed
            clip_ratio * image.std()
        double_check : bool, optional
            Optionally makes sure that the data contains only valid entries.
            Used for debugging.
        num_threads : int, optional
            If supplied, will use 'num_threads' number of OpenMP threads during
            the rendering.  Defaults to 0, which uses the environment variable
            OMP_NUM_THREADS.

        Returns
        -------
        image : array
            An (N,M,3) array of the final returned values, in float64 form.
        """
        if num_threads is None:
            num_threads=get_num_threads()
        image = self.new_image()
        args = self.get_sampler_args(image)
        sampler = self.get_sampler(args)
        self.initialize_source()
        image = ImageArray(self._render(double_check, num_threads, 
                                        image, sampler),
                           info=self.get_information())
        self.save_image(fn, clip_ratio, image)
        return image

    def show(self, clip_ratio = None):
        r"""This will take a snapshot and display the resultant image in the
        IPython notebook.

        If yt is being run from within an IPython session, and it is able to
        determine this, this function will snapshot and send the resultant
        image to the IPython notebook for display.

        If yt can't determine if it's inside an IPython session, it will raise
        YTNotInsideNotebook.

        Parameters
        ----------
        clip_ratio : float, optional
            If supplied, the 'max_val' argument to write_bitmap will be handed
            clip_ratio * image.std()

        Examples
        --------

        >>> cam.show()

        """
        if "__IPYTHON__" in dir(__builtin__):
            from IPython.core.displaypub import publish_display_data
            image = self.snapshot()
            if clip_ratio is not None: clip_ratio *= image.std()
            data = write_bitmap(image, None, clip_ratio)
            publish_display_data(
                'yt.visualization.volume_rendering.camera.Camera',
                {'image/png' : data}
            )
        else:
            raise YTNotInsideNotebook


    def set_default_light_dir(self):
        self.light_dir = [1.,1.,1.]

    def set_default_light_rgba(self):
        self.light_rgba = [1.,1.,1.,1.]

    def zoom(self, factor):
        r"""Change the distance to the focal point.

        This will zoom the camera in by some `factor` toward the focal point,
        along the current view direction, modifying the left/right and up/down
        extents as well.

        Parameters
        ----------
        factor : float
            The factor by which to reduce the distance to the focal point.


        Notes
        -----

        You will need to call snapshot() again to get a new image.

        """
        self.width = [w / factor for w in self.width]
        self._setup_box_properties(self.width, self.center, self.orienter.unit_vectors)

    def zoomin(self, final, n_steps, clip_ratio = None):
        r"""Loop over a zoomin and return snapshots along the way.

        This will yield `n_steps` snapshots until the current view has been
        zooming in to a final factor of `final`.

        Parameters
        ----------
        final : float
            The zoom factor, with respect to current, desired at the end of the
            sequence.
        n_steps : int
            The number of zoom snapshots to make.
        clip_ratio : float, optional
            If supplied, the 'max_val' argument to write_bitmap will be handed
            clip_ratio * image.std()


        Examples
        --------

        >>> for i, snapshot in enumerate(cam.zoomin(100.0, 10)):
        ...     iw.write_bitmap(snapshot, "zoom_%04i.png" % i)
        """
        f = final**(1.0/n_steps)
        for i in xrange(n_steps):
            self.zoom(f)
            yield self.snapshot(clip_ratio = clip_ratio)

    def move_to(self, final, n_steps, final_width=None, exponential=False, clip_ratio = None):
        r"""Loop over a look_at

        This will yield `n_steps` snapshots until the current view has been
        moved to a final center of `final` with a final width of final_width.

        Parameters
        ----------
        final : array_like
            The final center to move to after `n_steps`
        n_steps : int
            The number of look_at snapshots to make.
        final_width: float or array_like, optional
            Specifies the final width after `n_steps`.  Useful for
            moving and zooming at the same time.
        exponential : boolean
            Specifies whether the move/zoom transition follows an
            exponential path toward the destination or linear
        clip_ratio : float, optional
            If supplied, the 'max_val' argument to write_bitmap will be handed
            clip_ratio * image.std()
            
        Examples
        --------

        >>> for i, snapshot in enumerate(cam.move_to([0.2,0.3,0.6], 10)):
        ...     iw.write_bitmap(snapshot, "move_%04i.png" % i)
        """
        self.center = np.array(self.center)
        dW = None
        if exponential:
            if final_width is not None:
                if not iterable(final_width):
                    width = np.array([final_width, final_width, final_width]) 
                    # left/right, top/bottom, front/back 
                if (self.center == 0.0).all():
                    self.center += (np.array(final) - self.center) / (10. * n_steps)
                final_zoom = final_width/np.array(self.width)
                dW = final_zoom**(1.0/n_steps)
            else:
                dW = np.array([1.0,1.0,1.0])
            position_diff = (np.array(final)/self.center)*1.0
            dx = position_diff**(1.0/n_steps)
        else:
            if final_width is not None:
                if not iterable(final_width):
                    width = np.array([final_width, final_width, final_width]) 
                    # left/right, top/bottom, front/back
                dW = (1.0*final_width-np.array(self.width))/n_steps
            else:
                dW = np.array([0.0,0.0,0.0])
            dx = (np.array(final)-self.center)*1.0/n_steps
        for i in xrange(n_steps):
            if exponential:
                self.switch_view(center=self.center*dx, width=self.width*dW)
            else:
                self.switch_view(center=self.center+dx, width=self.width+dW)
            yield self.snapshot(clip_ratio = clip_ratio)

    def rotate(self, theta, rot_vector=None):
        r"""Rotate by a given angle

        Rotate the view.  If `rot_vector` is None, rotation will occur
        around the `north_vector`.

        Parameters
        ----------
        theta : float, in radians
             Angle (in radians) by which to rotate the view.
        rot_vector  : array_like, optional
            Specify the rotation vector around which rotation will
            occur.  Defaults to None, which sets rotation around
            `north_vector`

        Examples
        --------

        >>> cam.rotate(np.pi/4)
        """
        if rot_vector is None:
            rot_vector = self.rotation_vector
          
        R = get_rotation_matrix(theta, rot_vector)

        normal_vector = self.front_center-self.center

        self.switch_view(normal_vector=np.dot(R,normal_vector))

    def roll(self, theta):
        r"""Roll by a given angle

        Roll the view.

        Parameters
        ----------
        theta : float, in radians
             Angle (in radians) by which to roll the view.

        Examples
        --------

        >>> cam.roll(np.pi/4)
        """
        rot_vector = self.orienter.normal_vector
        R = get_rotation_matrix(theta, rot_vector)
        north_vector = self.orienter.unit_vectors[1]
        self.switch_view(north_vector=np.dot(R, north_vector))

    def rotation(self, theta, n_steps, rot_vector=None, clip_ratio = None):
        r"""Loop over rotate, creating a rotation

        This will yield `n_steps` snapshots until the current view has been
        rotated by an angle `theta`

        Parameters
        ----------
        theta : float, in radians
            Angle (in radians) by which to rotate the view.
        n_steps : int
            The number of look_at snapshots to make.
        rot_vector  : array_like, optional
            Specify the rotation vector around which rotation will
            occur.  Defaults to None, which sets rotation around the
            original `north_vector`
        clip_ratio : float, optional
            If supplied, the 'max_val' argument to write_bitmap will be handed
            clip_ratio * image.std()

        Examples
        --------

        >>> for i, snapshot in enumerate(cam.rotation(np.pi, 10)):
        ...     iw.write_bitmap(snapshot, 'rotation_%04i.png' % i)
        """

        dtheta = (1.0*theta)/n_steps
        for i in xrange(n_steps):
            self.rotate(dtheta, rot_vector=rot_vector)
            yield self.snapshot(clip_ratio = clip_ratio)

data_object_registry["camera"] = Camera

class InteractiveCamera(Camera):
    frames = []

    def snapshot(self, fn=None, clip_ratio=None):
        import matplotlib.pylab as pylab
        pylab.figure(2)
        self.transfer_function.show()
        pylab.draw()
        im = Camera.snapshot(self, fn, clip_ratio)
        pylab.figure(1)
        pylab.imshow(im / im.max())
        pylab.draw()
        self.frames.append(im)

    def rotation(self, theta, n_steps, rot_vector=None):
        for frame in Camera.rotation(self, theta, n_steps, rot_vector):
            if frame is not None:
                self.frames.append(frame)
                
    def zoomin(self, final, n_steps):
        for frame in Camera.zoomin(self, final, n_steps):
            if frame is not None:
                self.frames.append(frame)
                
    def clear_frames(self):
        del self.frames
        self.frames = []
        
    def save_frames(self, basename, clip_ratio=None):
        for i, frame in enumerate(self.frames):
            fn = basename + '_%04i.png'%i
            if clip_ratio is not None:
                write_bitmap(frame, fn, clip_ratio*image.std())
            else:
                write_bitmap(frame, fn)

data_object_registry["interactive_camera"] = InteractiveCamera

class PerspectiveCamera(Camera):
    expand_factor = 1.0
    def __init__(self, *args, **kwargs):
        self.expand_factor = kwargs.pop('expand_factor', 1.0)
        Camera.__init__(self, *args, **kwargs)

    def get_sampler_args(self, image):
        # We should move away from pre-generation of vectors like this and into
        # the usage of on-the-fly generation in the VolumeIntegrator module
        # We might have a different width and back_center
        dl = (self.back_center - self.front_center)
        self.front_center += self.expand_factor*dl
        self.back_center -= dl

        px = np.linspace(-self.width[0]/2.0, self.width[0]/2.0,
                         self.resolution[0])[:,None]
        py = np.linspace(-self.width[1]/2.0, self.width[1]/2.0,
                         self.resolution[1])[None,:]
        inv_mat = self.orienter.inv_mat
        positions = np.zeros((self.resolution[0], self.resolution[1], 3),
                          dtype='float64', order='C')
        positions[:,:,0] = inv_mat[0,0]*px+inv_mat[0,1]*py+self.back_center[0]
        positions[:,:,1] = inv_mat[1,0]*px+inv_mat[1,1]*py+self.back_center[1]
        positions[:,:,2] = inv_mat[2,0]*px+inv_mat[2,1]*py+self.back_center[2]
        bounds = (px.min(), px.max(), py.min(), py.max())

        # We are likely adding on an odd cutting condition here
        vectors = self.front_center - positions
        positions = self.front_center - 1.0*(((self.back_center-self.front_center)**2).sum())**0.5*vectors
        vectors = (self.front_center - positions)

        uv = np.ones(3, dtype='float64')
        image.shape = (self.resolution[0]**2,1,3)
        vectors.shape = (self.resolution[0]**2,1,3)
        positions.shape = (self.resolution[0]**2,1,3)
        args = (positions, vectors, self.back_center, 
                (0.0,1.0,0.0,1.0),
                image, uv, uv,
                np.zeros(3, dtype='float64'), 
                self.transfer_function, self.sub_samples)
        return args

    def _render(self, double_check, num_threads, image, sampler):
        pbar = get_pbar("Ray casting", (self.volume.brick_dimensions + 1).prod(axis=-1).sum())
        total_cells = 0
        if double_check:
            for brick in self.volume.bricks:
                for data in brick.my_data:
                    if np.any(np.isnan(data)):
                        raise RuntimeError

        view_pos = self.front_center
        for brick in self.volume.traverse(view_pos, self.front_center, image):
            sampler(brick, num_threads=num_threads)
            total_cells += np.prod(brick.my_data[0].shape)
            pbar.update(total_cells)

        pbar.finish()
        image = sampler.aimage
        self.finalize_image(image)
        return image


    def finalize_image(self, image):
        image.shape = self.resolution[0], self.resolution[0], 3

def corners(left_edge, right_edge):
    return np.array([
      [left_edge[:,0], left_edge[:,1], left_edge[:,2]],
      [right_edge[:,0], left_edge[:,1], left_edge[:,2]],
      [right_edge[:,0], right_edge[:,1], left_edge[:,2]],
      [right_edge[:,0], right_edge[:,1], right_edge[:,2]],
      [left_edge[:,0], right_edge[:,1], right_edge[:,2]],
      [left_edge[:,0], left_edge[:,1], right_edge[:,2]],
      [right_edge[:,0], left_edge[:,1], right_edge[:,2]],
      [left_edge[:,0], right_edge[:,1], left_edge[:,2]],
    ], dtype='float64')

class HEALpixCamera(Camera):

    _sampler_object = None 
    
    def __init__(self, center, radius, nside,
                 transfer_function = None, fields = None,
                 sub_samples = 5, log_fields = None, volume = None,
                 pf = None, use_kd=True, no_ghost=False, use_light=False):
        ParallelAnalysisInterface.__init__(self)
        if pf is not None: self.pf = pf
        self.center = np.array(center, dtype='float64')
        self.radius = radius
        self.nside = nside
        self.use_kd = use_kd
        if transfer_function is None:
            transfer_function = ProjectionTransferFunction()
        self.transfer_function = transfer_function

        if isinstance(self.transfer_function, ProjectionTransferFunction):
            self._sampler_object = ProjectionSampler
        else:
            self._sampler_object = VolumeRenderSampler

        if fields is None: fields = ["Density"]
        self.fields = fields
        self.sub_samples = sub_samples
        self.log_fields = log_fields
        self.use_light = use_light
        self.light_dir = None
        self.light_rgba = None
        if volume is None:
            volume = AMRKDTree(self.pf, fields=self.fields, no_ghost=no_ghost,
                               log_fields=log_fields)
        self.use_kd = isinstance(volume, AMRKDTree)
        self.volume = volume

    def new_image(self):
        image = np.zeros((12 * self.nside ** 2, 1, 3), dtype='float64', order='C')
        return image

    def get_sampler_args(self, image):
        nv = 12 * self.nside ** 2
        vs = arr_pix2vec_nest(self.nside, np.arange(nv))
        vs *= self.radius
        vs.shape = nv, 1, 3
        uv = np.ones(3, dtype='float64')
        positions = np.ones((nv, 1, 3), dtype='float64') * self.center
        args = (positions, vs, self.center,
                (0.0, 1.0, 0.0, 1.0),
                image, uv, uv,
                np.zeros(3, dtype='float64'),
                self.transfer_function, self.sub_samples)
        return args
 

    def _render(self, double_check, num_threads, image, sampler):
        pbar = get_pbar("Ray casting", (self.volume.brick_dimensions + 1).prod(axis=-1).sum())
        total_cells = 0
        if double_check:
            for brick in self.volume.bricks:
                for data in brick.my_data:
                    if np.any(np.isnan(data)):
                        raise RuntimeError
        
        view_pos = self.center
        for brick in self.volume.traverse(view_pos, None, image):
            sampler(brick, num_threads=num_threads)
            total_cells += np.prod(brick.my_data[0].shape)
            pbar.update(total_cells)
        
        pbar.finish()
        image = sampler.aimage

        self.finalize_image(image)

        return image

    def get_information(self):
        info_dict = {'fields':self.fields,
                     'type':self.__class__.__name__,
                     'center':self.center,
                     'radius':self.radius,
                     'dataset':self.pf.fullpath}
        return info_dict


    def snapshot(self, fn = None, clip_ratio = None, double_check = False,
                 num_threads = 0, clim = None, label = None):
        r"""Ray-cast the camera.

        This method instructs the camera to take a snapshot -- i.e., call the ray
        caster -- based on its current settings.

        Parameters
        ----------
        fn : string, optional
            If supplied, the image will be saved out to this before being
            returned.  Scaling will be to the maximum value.
        clip_ratio : float, optional
            If supplied, the 'max_val' argument to write_bitmap will be handed
            clip_ratio * image.std()

        Returns
        -------
        image : array
            An (N,M,3) array of the final returned values, in float64 form.
        """
        if num_threads is None:
            num_threads=get_num_threads()
        image = self.new_image()
        args = self.get_sampler_args(image)
        sampler = self.get_sampler(args)
        self.volume.initialize_source()
        image = ImageArray(self._render(double_check, num_threads, 
                                        image, sampler),
                           info=self.get_information())
        self.save_image(fn, clim, image, label = label)
        return image

    def save_image(self, fn, clim, image, label = None):
        if self.comm.rank is 0 and fn is not None:
            # This assumes Density; this is a relatively safe assumption.
            import matplotlib.figure
            import matplotlib.backends.backend_agg
            phi, theta = np.mgrid[0.0:2*np.pi:800j, 0:np.pi:800j]
            pixi = arr_ang2pix_nest(self.nside, theta.ravel(), phi.ravel())
            image *= self.radius * self.pf['cm']
            img = np.log10(image[:,0,0][pixi]).reshape((800,800))

            fig = matplotlib.figure.Figure((10, 5))
            ax = fig.add_subplot(1,1,1,projection='hammer')
            implot = ax.imshow(img, extent=(-np.pi,np.pi,-np.pi/2,np.pi/2), clip_on=False, aspect=0.5)
            cb = fig.colorbar(implot, orientation='horizontal')

            if label == None:
                cb.set_label("Projected %s" % self.fields[0])
            else:
                cb.set_label(label)
            if clim is not None: cb.set_clim(*clim)
            ax.xaxis.set_ticks(())
            ax.yaxis.set_ticks(())
            canvas = matplotlib.backends.backend_agg.FigureCanvasAgg(fig)
            canvas.print_figure(fn)


class AdaptiveHEALpixCamera(Camera):
    def __init__(self, center, radius, nside,
                 transfer_function = None, fields = None,
                 sub_samples = 5, log_fields = None, volume = None,
                 pf = None, use_kd=True, no_ghost=False,
                 rays_per_cell = 0.1, max_nside = 8192):
        ParallelAnalysisInterface.__init__(self)
        if pf is not None: self.pf = pf
        self.center = np.array(center, dtype='float64')
        self.radius = radius
        self.use_kd = use_kd
        if transfer_function is None:
            transfer_function = ProjectionTransferFunction()
        self.transfer_function = transfer_function
        if fields is None: fields = ["Density"]
        self.fields = fields
        self.sub_samples = sub_samples
        self.log_fields = log_fields
        if volume is None:
            volume = AMRKDTree(self.pf, fields=self.fields, no_ghost=no_ghost,
                               log_fields=log_fields)
        self.use_kd = isinstance(volume, AMRKDTree)
        self.volume = volume
        self.initial_nside = nside
        self.rays_per_cell = rays_per_cell
        self.max_nside = max_nside

    def snapshot(self, fn = None):
        tfp = TransferFunctionProxy(self.transfer_function)
        tfp.ns = self.sub_samples
        self.volume.initialize_source()
        mylog.info("Adaptively rendering.")
        pbar = get_pbar("Ray casting",
                        (self.volume.brick_dimensions + 1).prod(axis=-1).sum())
        total_cells = 0
        bricks = [b for b in self.volume.traverse(None, self.center, None)][::-1]
        left_edges = np.array([b.LeftEdge for b in bricks])
        right_edges = np.array([b.RightEdge for b in bricks])
        min_dx = min(((b.RightEdge[0] - b.LeftEdge[0])/b.my_data[0].shape[0]
                     for b in bricks))
        # We jitter a bit if we're on a boundary of our initial grid
        for i in range(3):
            if bricks[0].LeftEdge[i] == self.center[i]:
                self.center += 1e-2 * min_dx
            elif bricks[0].RightEdge[i] == self.center[i]:
                self.center -= 1e-2 * min_dx
        ray_source = AdaptiveRaySource(self.center, self.rays_per_cell,
                                       self.initial_nside, self.radius,
                                       bricks, left_edges, right_edges, self.max_nside)
        for i,brick in enumerate(bricks):
            ray_source.integrate_brick(brick, tfp, i, left_edges, right_edges,
                                       bricks)
            total_cells += np.prod(brick.my_data[0].shape)
            pbar.update(total_cells)
        pbar.finish()
        info, values = ray_source.get_rays()
        return info, values

class StereoPairCamera(Camera):
    def __init__(self, original_camera, relative_separation = 0.005):
        ParallelAnalysisInterface.__init__(self)
        self.original_camera = original_camera
        self.relative_separation = relative_separation

    def split(self):
        oc = self.original_camera
        uv = oc.orienter.unit_vectors
        c = oc.center
        fc = oc.front_center
        wx, wy, wz = oc.width
        left_normal = fc + uv[1] * 0.5*self.relative_separation * wx - c
        right_normal = fc - uv[1] * 0.5*self.relative_separation * wx - c
        left_camera = Camera(c, left_normal, oc.width,
                             oc.resolution, oc.transfer_function, north_vector=uv[0],
                             volume=oc.volume, fields=oc.fields, log_fields=oc.log_fields,
                             sub_samples=oc.sub_samples, pf=oc.pf)
        right_camera = Camera(c, right_normal, oc.width,
                             oc.resolution, oc.transfer_function, north_vector=uv[0],
                             volume=oc.volume, fields=oc.fields, log_fields=oc.log_fields,
                             sub_samples=oc.sub_samples, pf=oc.pf)
        return (left_camera, right_camera)

class FisheyeCamera(Camera):
    def __init__(self, center, radius, fov, resolution,
                 transfer_function = None, fields = None,
                 sub_samples = 5, log_fields = None, volume = None,
                 pf = None, no_ghost=False, rotation = None, use_light=False):
        ParallelAnalysisInterface.__init__(self)
        self.use_light = use_light
        self.light_dir = None
        self.light_rgba = None
        if rotation is None: rotation = np.eye(3)
        self.rotation_matrix = rotation
        if pf is not None: self.pf = pf
        self.center = np.array(center, dtype='float64')
        self.radius = radius
        self.fov = fov
        if iterable(resolution):
            raise RuntimeError("Resolution must be a single int")
        self.resolution = resolution
        if transfer_function is None:
            transfer_function = ProjectionTransferFunction()
        self.transfer_function = transfer_function
        if fields is None: fields = ["Density"]
        self.fields = fields
        self.sub_samples = sub_samples
        self.log_fields = log_fields
        if volume is None:
            volume = AMRKDTree(self.pf, fields=self.fields, no_ghost=no_ghost,
                               log_fields=log_fields)
        self.volume = volume

    def new_image(self):
        image = np.zeros((self.resolution**2,1,3), dtype='float64', order='C')
        return image
        
    def get_sampler_args(self, image):
        vp = arr_fisheye_vectors(self.resolution, self.fov)
        vp.shape = (self.resolution**2,1,3)
        vp2 = vp.copy()
        for i in range(3):
            vp[:,:,i] = (vp2 * self.rotation_matrix[:,i]).sum(axis=2)
        del vp2
        vp *= self.radius
        uv = np.ones(3, dtype='float64')
        positions = np.ones((self.resolution**2, 1, 3), dtype='float64') * self.center

        args = (positions, vp, self.center,
                (0.0, 1.0, 0.0, 1.0),
                image, uv, uv,
                np.zeros(3, dtype='float64'),
                self.transfer_function, self.sub_samples)
        return args


    def finalize_image(self, image):
        image.shape = self.resolution, self.resolution, 3

    def _render(self, double_check, num_threads, image, sampler):
        pbar = get_pbar("Ray casting", (self.volume.brick_dimensions + 1).prod(axis=-1).sum())
        total_cells = 0
        if double_check:
            for brick in self.volume.bricks:
                for data in brick.my_data:
                    if np.any(np.isnan(data)):
                        raise RuntimeError
        
        view_pos = self.center
        for brick in self.volume.traverse(view_pos, None, image):
            sampler(brick, num_threads=num_threads)
            total_cells += np.prod(brick.my_data[0].shape)
            pbar.update(total_cells)
        
        pbar.finish()
        image = sampler.aimage

        self.finalize_image(image)

        return image

class MosaicFisheyeCamera(Camera):
    def __init__(self, center, radius, fov, resolution, focal_center=None,
                 transfer_function=None, fields=None,
                 sub_samples=5, log_fields=None, volume=None,
                 pf=None, l_max=None, no_ghost=False,nimx=1, nimy=1, procs_per_wg=None,
                 rotation=None):
        r"""A fisheye lens camera, taking adantage of image plane decomposition
        for parallelism..

        The camera represents the eye of an observer, which will be used to
        generate ray-cast volume renderings of the domain. In this case, the
        rays are defined by a fisheye lens

        Parameters
        ----------
        center : array_like
            The current "center" of the observer, from which the rays will be
            cast
        radius : float
            The radial distance to cast to
        resolution : int
            The number of pixels in each direction.  Must be a single int.
        volume : `yt.extensions.volume_rendering.HomogenizedVolume`, optional
            The volume to ray cast through.  Can be specified for finer-grained
            control, but otherwise will be automatically generated.
        fields : list of fields, optional
            This is the list of fields we want to volume render; defaults to
            Density.
        log_fields : list of bool, optional
            Whether we should take the log of the fields before supplying them to
            the volume rendering mechanism.
        sub_samples : int, optional
            The number of samples to take inside every cell per ray.
        pf : `~yt.data_objects.api.StaticOutput`
            For now, this is a require parameter!  But in the future it will become
            optional.  This is the parameter file to volume render.
        l_max: int, optional
            Specifies the maximum level to be rendered.  Also
            specifies the maximum level used in the AMRKDTree
            construction.  Defaults to None (all levels), and only
            applies if use_kd=True.
        no_ghost: bool, optional
            Optimization option.  If True, homogenized bricks will
            extrapolate out from grid instead of interpolating from
            ghost zones that have to first be calculated.  This can
            lead to large speed improvements, but at a loss of
            accuracy/smoothness in resulting image.  The effects are
            less notable when the transfer function is smooth and
            broad. Default: False
        nimx: int, optional
            The number by which to decompose the image plane into in the x
            direction.  Must evenly divide the resolution.
        nimy: int, optional
            The number by which to decompose the image plane into in the y 
            direction.  Must evenly divide the resolution.
        procs_per_wg: int, optional
            The number of processors to use on each sub-image. Within each
            subplane, the volume will be decomposed using the AMRKDTree with
            procs_per_wg processors.  

        Notes
        -----
            The product of nimx*nimy*procs_per_wg must be equal to or less than
            the total number of mpi processes.  

            Unlike the non-Mosaic camera, this will only return each sub-image
            to the root processor of each sub-image workgroup in order to save
            memory.  To save the final image, one must then call
            MosaicFisheyeCamera.save_image('filename')

        Examples
        --------

        >>> from yt.mods import *
        
        >>> pf = load('DD1717')
        
        >>> N = 512 # Pixels (1024^2)
        >>> c = (pf.domain_right_edge + pf.domain_left_edge)/2. # Center
        >>> radius = (pf.domain_right_edge - pf.domain_left_edge)/2.
        >>> fov = 180.0
        
        >>> field='Density'
        >>> mi,ma = pf.h.all_data().quantities['Extrema']('Density')[0]
        >>> mi,ma = np.log10(mi), np.log10(ma)
        
        # You may want to comment out the above lines and manually set the min and max
        # of the log of the Density field. For example:
        # mi,ma = -30.5,-26.5
        
        # Another good place to center the camera is close to the maximum density.
        # v,c = pf.h.find_max('Density')
        # c -= 0.1*radius
        
       
        # Construct transfer function
        >>> tf = ColorTransferFunction((mi-1, ma+1),nbins=1024)
        
        # Sample transfer function with Nc gaussians.  Use col_bounds keyword to limit
        # the color range to the min and max values, rather than the transfer function
        # bounds.
        >>> Nc = 5
        >>> tf.add_layers(Nc,w=0.005, col_bounds = (mi,ma), alpha=np.logspace(-2,0,Nc),
        >>>         colormap='RdBu_r')
        >>> 
        # Create the camera object. Use the keyword: no_ghost=True if a lot of time is
        # spent creating vertex-centered data. In this case I'm running with 8
        # processors, and am splitting the image plane into 4 pieces and using 2
        # processors on each piece.
        >>> cam = MosaicFisheyeCamera(c, radius, fov, N,
        >>>         transfer_function = tf, 
        >>>         sub_samples = 5, 
        >>>         pf=pf, 
        >>>         nimx=2,nimy=2,procs_per_wg=2)
        
        # Take a snapshot
        >>> im = cam.snapshot()
        
        # Save the image
        >>> cam.save_image('fisheye_mosaic.png')

        """

        ParallelAnalysisInterface.__init__(self)
        self.image_decomp = self.comm.size>1
        if self.image_decomp:
            PP = ProcessorPool()
            npatches = nimy*nimx
            if procs_per_wg is None:
                if (PP.size % npatches):
                    raise RuntimeError("Cannot evenly divide %i procs to %i patches" % (PP.size,npatches))
                else:
                    procs_per_wg = PP.size / npatches
            if (PP.size != npatches*procs_per_wg):
               raise RuntimeError("You need %i processors to utilize %i procs per one patch in [%i,%i] grid" 
                     % (npatches*procs_per_wg,procs_per_wg,nimx,nimy))
 
            for j in range(nimy):
                for i in range(nimx):
                    PP.add_workgroup(size=procs_per_wg, name='%04i_%04i'%(i,j))
                    
            for wg in PP.workgroups:
                if self.comm.rank in wg.ranks:
                    my_wg = wg
            
            self.global_comm = self.comm
            self.comm = my_wg.comm
            self.wg = my_wg
            self.imi = int(self.wg.name[0:4])
            self.imj = int(self.wg.name[5:9])
            mylog.info('My new communicator has the name %s' % self.wg.name)
            self.nimx = nimx
            self.nimy = nimy
        else:
            self.imi = 0
            self.imj = 0
            self.nimx = 1
            self.nimy = 1
        if pf is not None: self.pf = pf
        
        if rotation is None: rotation = np.eye(3)
        self.rotation_matrix = rotation
        
        self.normal_vector = np.array([0.,0.,1])
        self.north_vector = np.array([1.,0.,0.])
        self.east_vector = np.array([0.,1.,0.])
        self.rotation_vector = self.north_vector

        if iterable(resolution):
            raise RuntimeError("Resolution must be a single int")
        self.resolution = resolution
        self.center = np.array(center, dtype='float64')
        self.focal_center = focal_center
        self.radius = radius
        self.fov = fov
        if transfer_function is None:
            transfer_function = ProjectionTransferFunction()
        self.transfer_function = transfer_function
        if fields is None: fields = ["Density"]
        self.fields = fields
        self.sub_samples = sub_samples
        self.log_fields = log_fields
        if volume is None:
            volume = AMRKDTree(self.pf, fields=self.fields, no_ghost=no_ghost,
                               log_fields=log_fields,l_max=l_max)
        self.volume = volume
        self.vp = None
        self.image = None 

    def get_vector_plane(self):
        if self.focal_center is not None:
            rvec =  np.array(self.focal_center) - np.array(self.center)
            rvec /= (rvec**2).sum()**0.5
            angle = np.arccos( (self.normal_vector*rvec).sum()/( (self.normal_vector**2).sum()**0.5 *
                (rvec**2).sum()**0.5))
            rot_vector = np.cross(rvec, self.normal_vector)
            rot_vector /= (rot_vector**2).sum()**0.5
            
            self.rotation_matrix = get_rotation_matrix(angle,rot_vector)
            self.normal_vector = np.dot(self.rotation_matrix,self.normal_vector)
            self.north_vector = np.dot(self.rotation_matrix,self.north_vector)
            self.east_vector = np.dot(self.rotation_matrix,self.east_vector)
        else:
            self.focal_center = self.center + self.radius*self.normal_vector  
        dist = ((self.focal_center - self.center)**2).sum()**0.5
        # We now follow figures 4-7 of:
        # http://paulbourke.net/miscellaneous/domefisheye/fisheye/
        # ...but all in Cython.
        
        self.vp = arr_fisheye_vectors(self.resolution, self.fov, self.nimx, 
                self.nimy, self.imi, self.imj)
        
        self.vp = rotate_vectors(self.vp, self.rotation_matrix)

        self.center = self.focal_center - dist*self.normal_vector
        self.vp *= self.radius
        nx, ny = self.vp.shape[0], self.vp.shape[1]
        self.vp.shape = (nx*ny,1,3)

    def snapshot(self):
        if self.vp is None:
            self.get_vector_plane()

        nx,ny = self.resolution/self.nimx, self.resolution/self.nimy
        image = np.zeros((nx*ny,1,3), dtype='float64', order='C')
        uv = np.ones(3, dtype='float64')
        positions = np.ones((nx*ny, 1, 3), dtype='float64') * self.center
        vector_plane = VectorPlane(positions, self.vp, self.center,
                        (0.0, 1.0, 0.0, 1.0), image, uv, uv)
        tfp = TransferFunctionProxy(self.transfer_function)
        tfp.ns = self.sub_samples
        self.volume.initialize_source()
        mylog.info("Rendering fisheye of %s^2", self.resolution)
        pbar = get_pbar("Ray casting",
                        (self.volume.brick_dimensions + 1).prod(axis=-1).sum())

        total_cells = 0
        for brick in self.volume.traverse(None, self.center, image):
            brick.cast_plane(tfp, vector_plane)
            total_cells += np.prod(brick.my_data[0].shape)
            pbar.update(total_cells)
        pbar.finish()
        image.shape = (nx, ny, 3)

        if self.image is not None:
            del self.image
        image = ImageArray(image,
                           info=self.get_information())
        self.image = image
        return image

    def save_image(self, fn, clip_ratio=None):
        if '.png' not in fn:
            fn = fn + '.png'
        
        try:
            image = self.image
        except:
            mylog.error('You must first take a snapshot')
            raise(UserWarning)
        
        image = self.image
        nx,ny = self.resolution/self.nimx, self.resolution/self.nimy
        if self.image_decomp:
            if self.comm.rank == 0:
                if self.global_comm.rank == 0:
                    final_image = np.empty((nx*self.nimx, 
                        ny*self.nimy, 3),
                        dtype='float64',order='C')
                    final_image[:nx, :ny, :] = image
                    for j in range(self.nimy):
                        for i in range(self.nimx):
                            if i==0 and j==0: continue
                            arr = self.global_comm.recv_array((self.wg.size)*(j*self.nimx + i), tag = (self.wg.size)*(j*self.nimx + i))

                            final_image[i*nx:(i+1)*nx, j*ny:(j+1)*ny,:] = arr
                            del arr
                    if clip_ratio is not None:
                        write_bitmap(final_image, fn, clip_ratio*final_image.std())
                    else:
                        write_bitmap(final_image, fn)
                else:
                    self.global_comm.send_array(image, 0, tag = self.global_comm.rank)
        else:
            if self.comm.rank == 0:
                if clip_ratio is not None:
                    write_bitmap(image, fn, clip_ratio*image.std())
                else:
                    write_bitmap(image, fn)
        return

    def rotate(self, theta, rot_vector=None, keep_focus=True):
        r"""Rotate by a given angle

        Rotate the view.  If `rot_vector` is None, rotation will occur
        around the `north_vector`.

        Parameters
        ----------
        theta : float, in radians
             Angle (in radians) by which to rotate the view.
        rot_vector  : array_like, optional
            Specify the rotation vector around which rotation will
            occur.  Defaults to None, which sets rotation around
            `north_vector`

        Examples
        --------

        >>> cam.rotate(np.pi/4)
        """
        if rot_vector is None:
            rot_vector = self.north_vector

        dist = ((self.focal_center - self.center)**2).sum()**0.5

        R = get_rotation_matrix(theta, rot_vector)

        self.vp = rotate_vectors(self.vp, R)
        self.normal_vector = np.dot(R,self.normal_vector)
        self.north_vector = np.dot(R,self.north_vector)
        self.east_vector = np.dot(R,self.east_vector)

        if keep_focus:
            self.center = self.focal_center - dist*self.normal_vector

    def rotation(self, theta, n_steps, rot_vector=None, keep_focus=True):
        r"""Loop over rotate, creating a rotation

        This will yield `n_steps` snapshots until the current view has been
        rotated by an angle `theta`

        Parameters
        ----------
        theta : float, in radians
            Angle (in radians) by which to rotate the view.
        n_steps : int
            The number of look_at snapshots to make.
        rot_vector  : array_like, optional
            Specify the rotation vector around which rotation will
            occur.  Defaults to None, which sets rotation around the
            original `north_vector`

        Examples
        --------

        >>> for i, snapshot in enumerate(cam.rotation(np.pi, 10)):
        ...     iw.write_bitmap(snapshot, 'rotation_%04i.png' % i)
        """

        dtheta = (1.0*theta)/n_steps
        for i in xrange(n_steps):
            self.rotate(dtheta, rot_vector=rot_vector, keep_focus=keep_focus)
            yield self.snapshot()

    def move_to(self,final,n_steps,exponential=False):
        r"""Loop over a look_at

        This will yield `n_steps` snapshots until the current view has been
        moved to a final center of `final`.

        Parameters
        ----------
        final : array_like
            The final center to move to after `n_steps`
        n_steps : int
            The number of look_at snapshots to make.
        exponential : boolean
            Specifies whether the move/zoom transition follows an
            exponential path toward the destination or linear

        Examples
        --------

        >>> for i, snapshot in enumerate(cam.move_to([0.2,0.3,0.6], 10)):
        ...     cam.save_image('move_%04i.png' % i)
        """
        if exponential:
            position_diff = (np.array(final)/self.center)*1.0
            dx = position_diff**(1.0/n_steps)
        else:
            dx = (np.array(final) - self.center)*1.0/n_steps
        for i in xrange(n_steps):
            if exponential:
                self.center *= dx
            else:
                self.center += dx
            yield self.snapshot()

def allsky_projection(pf, center, radius, nside, field, weight = None,
                      inner_radius = 10, rotation = None):
    r"""Project through a parameter file, through an allsky-method
    decomposition from HEALpix, and return the image plane.

    This function will accept the necessary items to integrate through a volume
    over 4pi and return the integrated field of view to the user.  Note that if
    a weight is supplied, it will multiply the pre-interpolated values
    together.

    Parameters
    ----------
    pf : `~yt.data_objects.api.StaticOutput`
        This is the parameter file to volume render.
    center : array_like
        The current "center" of the view port -- the focal point for the
        camera.
    radius : float or list of floats
        The radius to integrate out to of the image.
    nside : int
        The HEALpix degree.  The number of rays integrated is 12*(Nside**2)
        Must be a power of two!
    field : string
        The field to project through the volume
    weight : optional, default None
        If supplied, the field will be pre-multiplied by this, then divided by
        the integrated value of this field.  This returns an average rather
        than a sum.
    inner_radius : optional, float, defaults to 0.05
        The radius of the inner clipping plane, in units of the dx at the point
        at which the volume rendering is centered.  This avoids unphysical
        effects of nearby cells.
    rotation : optional, 3x3 array
        If supplied, the vectors will be rotated by this.  You can construct
        this by, for instance, calling np.array([v1,v2,v3]) where those are the
        three reference planes of an orthogonal frame (see ortho_find).

    Returns
    -------
    image : array
        An ((Nside**2)*12,1,3) array of the final integrated values, in float64 form.

    Examples
    --------

    >>> image = allsky_projection(pf, [0.5, 0.5, 0.5], 1.0/pf['mpc'],
                      32, "Temperature", "Density")
    >>> plot_allsky_healpix(image, 32, "healpix.png")

    """
    # We manually modify the ProjectionTransferFunction to get it to work the
    # way we want, with a second field that's also passed through.
    fields = [field]
    center = np.array(center, dtype='float64')
    if weight is not None:
        # This is a temporary field, which we will remove at the end.
        def _make_wf(f, w):
            def temp_weightfield(a, b):
                tr = b[f].astype("float64") * b[w]
                return tr
            return temp_weightfield
        pf.field_info.add_field("temp_weightfield",
            function=_make_wf(field, weight))
        fields = ["temp_weightfield", weight]
    nv = 12*nside**2
    image = np.zeros((nv,1,3), dtype='float64', order='C')
    vs = arr_pix2vec_nest(nside, np.arange(nv))
    vs.shape = (nv,1,3)
    if rotation is not None:
        vs2 = vs.copy()
        for i in range(3):
            vs[:,:,i] = (vs2 * rotation[:,i]).sum(axis=2)
    else:
        vs += 1e-8
    positions = np.ones((nv, 1, 3), dtype='float64', order='C') * center
    dx = min(g.dds.min() for g in pf.h.find_point(center)[0])
    positions += inner_radius * dx * vs
    vs *= radius
    uv = np.ones(3, dtype='float64')
    grids = pf.h.sphere(center, radius)._grids
    sampler = ProjectionSampler(positions, vs, center, (0.0, 0.0, 0.0, 0.0),
                                image, uv, uv, np.zeros(3, dtype='float64'))
    pb = get_pbar("Sampling ", len(grids))
    for i,grid in enumerate(grids):
        data = [grid[field] * grid.child_mask.astype('float64')
                for field in fields]
        pg = PartitionedGrid(
            grid.id, data,
            grid.LeftEdge, grid.RightEdge,
            grid.ActiveDimensions.astype("int64"))
        grid.clear_data()
        sampler(pg)
        pb.update(i)
    pb.finish()
    image = sampler.aimage
    if weight is None:
        dl = radius * pf.units[pf.field_info[field].projection_conversion]
        image *= dl
    else:
        image[:,:,0] /= image[:,:,1]
        pf.field_info.pop("temp_weightfield")
        for g in pf.h.grids:
            if "temp_weightfield" in g.keys():
                del g["temp_weightfield"]
    return image[:,0,0]

def plot_allsky_healpix(image, nside, fn, label = "", rotation = None,
                        take_log = True, resolution=512, cmin=None, cmax=None):
    import matplotlib.figure
    import matplotlib.backends.backend_agg
    if rotation is None: rotation = np.eye(3).astype("float64")

    img, count = pixelize_healpix(nside, image, resolution, resolution, rotation)

    fig = matplotlib.figure.Figure((10, 5))
    ax = fig.add_subplot(1,1,1,projection='aitoff')
    if take_log: func = np.log10
    else: func = lambda a: a
    implot = ax.imshow(func(img), extent=(-np.pi,np.pi,-np.pi/2,np.pi/2),
                       clip_on=False, aspect=0.5, vmin=cmin, vmax=cmax)
    cb = fig.colorbar(implot, orientation='horizontal')
    cb.set_label(label)
    ax.xaxis.set_ticks(())
    ax.yaxis.set_ticks(())
    canvas = matplotlib.backends.backend_agg.FigureCanvasAgg(fig)
    canvas.print_figure(fn)
    return img, count

class ProjectionCamera(Camera):
    def __init__(self, center, normal_vector, width, resolution,
            field, weight=None, volume=None, no_ghost = False, 
            le=None, re=None,
            north_vector=None, pf=None, interpolated=False):

        if not interpolated:
            volume = 1

        self.interpolated = interpolated
        self.field = field
        self.weight = weight
        self.resolution = resolution

        fields = [field]
        if self.weight is not None:
            # This is a temporary field, which we will remove at the end.
            def _make_wf(f, w):
                def temp_weightfield(a, b):
                    tr = b[f].astype("float64") * b[w]
                    return tr
                return temp_weightfield
            pf.field_info.add_field("temp_weightfield",
                function=_make_wf(self.field, self.weight))
            fields = ["temp_weightfield", self.weight]
        
        self.fields = fields
        self.log_fields = [False]*len(self.fields)
        Camera.__init__(self, center, normal_vector, width, resolution, None,
                fields = fields, pf=pf, volume=volume,
                log_fields=self.log_fields, 
                le=le, re=re, north_vector=north_vector,
                no_ghost=no_ghost)

    def get_sampler(self, args):
        if self.interpolated:
            sampler = InterpolatedProjectionSampler(*args)
        else:
            sampler = ProjectionSampler(*args)
        return sampler

    def initialize_source(self):
        if self.interpolated:
            Camera.initialize_source(self)
        else:
            pass

    def get_sampler_args(self, image):
        rotp = np.concatenate([self.orienter.inv_mat.ravel('F'), self.back_center.ravel()])
        args = (rotp, self.box_vectors[2], self.back_center,
            (-self.width[0]/2, self.width[0]/2,
             -self.width[1]/2, self.width[1]/2),
            image, self.orienter.unit_vectors[0], self.orienter.unit_vectors[1],
                np.array(self.width), self.sub_samples)
        return args

    def finalize_image(self,image):
        pf = self.pf
        if self.weight is None:
            dl = self.width[2] * pf.units[pf.field_info[self.field].projection_conversion]
            image *= dl
        else:
            image[:,:,0] /= image[:,:,1]
        return image[:,:,0]


    def _render(self, double_check, num_threads, image, sampler):
        # Calculate the eight corners of the box
        # Back corners ...
        if self.interpolated:
            return Camera._render(self, double_check, num_threads, image,
                    sampler)
        pf = self.pf
        width = self.width[2]
        north_vector = self.orienter.unit_vectors[0]
        east_vector = self.orienter.unit_vectors[1]
        normal_vector = self.orienter.unit_vectors[2]
        fields = self.fields

        mi = pf.domain_right_edge.copy()
        ma = pf.domain_left_edge.copy()
        for off1 in [-1, 1]:
            for off2 in [-1, 1]:
                for off3 in [-1, 1]:
                    this_point = (self.center + width/2. * off1 * north_vector
                                         + width/2. * off2 * east_vector
                                         + width/2. * off3 * normal_vector)
                    np.minimum(mi, this_point, mi)
                    np.maximum(ma, this_point, ma)
        # Now we have a bounding box.
        grids = pf.h.region(self.center, mi, ma)._grids

        pb = get_pbar("Sampling ", len(grids))
        for i,grid in enumerate(grids):
            data = [(grid[field] * grid.child_mask).astype("float64")
                    for field in fields]
            pg = PartitionedGrid(
                grid.id, data,
                grid.LeftEdge, grid.RightEdge, grid.ActiveDimensions.astype("int64"))
            grid.clear_data()
            sampler(pg, num_threads = num_threads)
            pb.update(i)
        pb.finish()

        image = sampler.aimage
        self.finalize_image(image)
        return image

    def save_image(self, fn, clip_ratio, image):
        if self.pf.field_info[self.field].take_log:
            im = np.log10(image)
        else:
            im = image
        if self.comm.rank is 0 and fn is not None:
            if clip_ratio is not None:
                write_image(im, fn)
            else:
                write_image(im, fn)

    def snapshot(self, fn = None, clip_ratio = None, double_check = False,
                 num_threads = 0):

        if num_threads is None:
            num_threads=get_num_threads()

        fields = [self.field]
        resolution = self.resolution

        image = self.new_image()

        args = self.get_sampler_args(image)

        sampler = self.get_sampler(args)

        self.initialize_source()

        image = ImageArray(self._render(double_check, num_threads, 
                                        image, sampler),
                           info=self.get_information())

        self.save_image(fn, clip_ratio, image)

        return image
    snapshot.__doc__ = Camera.snapshot.__doc__

data_object_registry["projection_camera"] = ProjectionCamera

def off_axis_projection(pf, center, normal_vector, width, resolution,
                        field, weight = None, 
                        volume = None, no_ghost = False, interpolated = False,
                        north_vector = None):
    r"""Project through a parameter file, off-axis, and return the image plane.

    This function will accept the necessary items to integrate through a volume
    at an arbitrary angle and return the integrated field of view to the user.
    Note that if a weight is supplied, it will multiply the pre-interpolated
    values together, then create cell-centered values, then interpolate within
    the cell to conduct the integration.

    Parameters
    ----------
    pf : `~yt.data_objects.api.StaticOutput`
        This is the parameter file to volume render.
    center : array_like
        The current 'center' of the view port -- the focal point for the
        camera.
    normal_vector : array_like
        The vector between the camera position and the center.
    width : float or list of floats
        The current width of the image.  If a single float, the volume is
        cubical, but if not, it is left/right, top/bottom, front/back
    resolution : int or list of ints
        The number of pixels in each direction.
    field : string
        The field to project through the volume
    weight : optional, default None
        If supplied, the field will be pre-multiplied by this, then divided by
        the integrated value of this field.  This returns an average rather
        than a sum.
    volume : `yt.extensions.volume_rendering.HomogenizedVolume`, optional
        The volume to ray cast through.  Can be specified for finer-grained
        control, but otherwise will be automatically generated.
    no_ghost: bool, optional
        Optimization option.  If True, homogenized bricks will
        extrapolate out from grid instead of interpolating from
        ghost zones that have to first be calculated.  This can
        lead to large speed improvements, but at a loss of
        accuracy/smoothness in resulting image.  The effects are
        less notable when the transfer function is smooth and
        broad. Default: True
    interpolated : optional, default False
        If True, the data is first interpolated to vertex-centered data, 
        then tri-linearly interpolated along the ray. Not suggested for 
        quantitative studies.

    Returns
    -------
    image : array
        An (N,N) array of the final integrated values, in float64 form.

    Examples
    --------

    >>> image = off_axis_projection(pf, [0.5, 0.5, 0.5], [0.2,0.3,0.4],
                      0.2, N, "Temperature", "Density")
    >>> write_image(np.log10(image), "offaxis.png")

    """
    projcam = ProjectionCamera(center, normal_vector, width, resolution,
                               field, weight=weight, pf=pf, volume=volume,
                               no_ghost=no_ghost, interpolated=interpolated, 
                               north_vector=north_vector)
    image = projcam.snapshot()
    if weight is not None:
        pf.field_info.pop("temp_weightfield")
    del projcam
    return image[:,:,0]