Source

yt-3.0 / yt / visualization / plot_modifications.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
""" 
Callbacks to add additional functionality on to plots.

Author: Matthew Turk <matthewturk@gmail.com>
Affiliation: KIPAC/SLAC/Stanford
Author: J. S. Oishi <jsoishi@astro.berkeley.edu>
Affiliation: UC Berkeley
Author: Stephen Skory <s@skory.us>
Affiliation: UC San Diego
Author: Anthony Scopatz <scopatz@gmail.com>
Affiliation: The University of Chicago
Homepage: http://yt-project.org/
License:
  Copyright (C) 2008-2012 Matthew Turk, JS Oishi, Stephen Skory, Anthony Scopatz.  
  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 numpy as np

from yt.funcs import *
from _mpl_imports import *
from yt.utilities.definitions import \
    x_dict, x_names, \
    y_dict, y_names, \
    axis_names, \
    axis_labels
import _MPL

callback_registry = {}

class PlotCallback(object):
    class __metaclass__(type):
        def __init__(cls, name, b, d):
            type.__init__(cls, name, b, d)
            callback_registry[name] = cls

    def __init__(self, *args, **kwargs):
        pass

    def convert_to_plot(self, plot, coord, offset = True):
        # coord should be a 2 x ncoord array-like datatype.
        try:
            ncoord = np.array(coord).shape[1]
        except IndexError:
            ncoord = 1

        # Convert the data and plot limits to tiled numpy arrays so that
        # convert_to_plot is automatically vectorized.

        x0 = np.tile(plot.xlim[0],ncoord)
        x1 = np.tile(plot.xlim[1],ncoord)
        xx0 = np.tile(plot._axes.get_xlim()[0],ncoord)
        xx1 = np.tile(plot._axes.get_xlim()[1],ncoord)
        
        y0 = np.tile(plot.ylim[0],ncoord)
        y1 = np.tile(plot.ylim[1],ncoord)
        yy0 = np.tile(plot._axes.get_ylim()[0],ncoord)
        yy1 = np.tile(plot._axes.get_ylim()[1],ncoord)
        
        # We need a special case for when we are only given one coordinate.
        if np.array(coord).shape == (2,):
            return ((coord[0]-x0)/(x1-x0)*(xx1-xx0) + xx0,
                    (coord[1]-y0)/(y1-y0)*(yy1-yy0) + yy0)
        else:
            return ((coord[0][:]-x0)/(x1-x0)*(xx1-xx0) + xx0,
                    (coord[1][:]-y0)/(y1-y0)*(yy1-yy0) + yy0)

    def pixel_scale(self,plot):
        x0, x1 = plot.xlim
        xx0, xx1 = plot._axes.get_xlim()
        dx = (xx0 - xx1)/(x1 - x0)
        
        y0, y1 = plot.ylim
        yy0, yy1 = plot._axes.get_ylim()
        dy = (yy0 - yy1)/(y1 - y0)

        return (dx,dy)


class VelocityCallback(PlotCallback):
    _type_name = "velocity"
    def __init__(self, factor=16, scale=None, scale_units=None, normalize=False):
        """
        annotate_velocity(factor=16, scale=None, scale_units=None, normalize=False):
        
        Adds a 'quiver' plot of velocity to the plot, skipping all but
        every *factor* datapoint. *scale* is the data units per arrow
        length unit using *scale_units* (see
        matplotlib.axes.Axes.quiver for more info). if *normalize* is
        True, the velocity fields will be scaled by their local
        (in-plane) length, allowing morphological features to be more
        clearly seen for fields with substantial variation in field
        strength (normalize is not implemented and thus ignored for
        Cutting Planes).
        """
        PlotCallback.__init__(self)
        self.factor = factor
        self.scale  = scale
        self.scale_units = scale_units
        self.normalize = normalize

    def __call__(self, plot):
        # Instantiation of these is cheap
        if plot._type_name == "CuttingPlane":
            qcb = CuttingQuiverCallback("CuttingPlaneVelocityX",
                                        "CuttingPlaneVelocityY",
                                        self.factor)
        else:
            xv = "%s-velocity" % (x_names[plot.data.axis])
            yv = "%s-velocity" % (y_names[plot.data.axis])
            qcb = QuiverCallback(xv, yv, self.factor, scale=self.scale, scale_units=self.scale_units, normalize=self.normalize)
        return qcb(plot)

class MagFieldCallback(PlotCallback):
    _type_name = "magnetic_field"
    def __init__(self, factor=16, scale=None, scale_units=None, normalize=False):
        """
        annotate_magnetic_field(factor=16, scale=None, scale_units=None, normalize=False):

        Adds a 'quiver' plot of magnetic field to the plot, skipping all but
        every *factor* datapoint. *scale* is the data units per arrow
        length unit using *scale_units* (see
        matplotlib.axes.Axes.quiver for more info). if *normalize* is
        True, the magnetic fields will be scaled by their local
        (in-plane) length, allowing morphological features to be more
        clearly seen for fields with substantial variation in field strength.
        """
        PlotCallback.__init__(self)
        self.factor = factor
        self.scale  = scale
        self.scale_units = scale_units
        self.normalize = normalize

    def __call__(self, plot):
        # Instantiation of these is cheap
        if plot._type_name == "CuttingPlane":
            qcb = CuttingQuiverCallback("CuttingPlaneBx",
                                        "CuttingPlaneBy",
                                        self.factor)
        else:
            xv = "B%s" % (x_names[plot.data.axis])
            yv = "B%s" % (y_names[plot.data.axis])
            qcb = QuiverCallback(xv, yv, self.factor, scale=self.scale, scale_units=self.scale_units, normalize=self.normalize)
        return qcb(plot)

class QuiverCallback(PlotCallback):
    _type_name = "quiver"
    def __init__(self, field_x, field_y, factor=16, scale=None, scale_units=None, normalize=False):
        """
        annotate_quiver(field_x, field_y, factor, scale=None, scale_units=None, normalize=False):

        Adds a 'quiver' plot to any plot, using the *field_x* and *field_y*
        from the associated data, skipping every *factor* datapoints
        *scale* is the data units per arrow length unit using *scale_units* 
        (see matplotlib.axes.Axes.quiver for more info)
        """
        PlotCallback.__init__(self)
        self.field_x = field_x
        self.field_y = field_y
        self.bv_x = self.bv_y = 0
        self.factor = factor
        self.scale = scale
        self.scale_units = scale_units
        self.normalize = normalize

    def __call__(self, plot):
        x0, x1 = plot.xlim
        y0, y1 = plot.ylim
        xx0, xx1 = plot._axes.get_xlim()
        yy0, yy1 = plot._axes.get_ylim()
        plot._axes.hold(True)
        nx = plot.image._A.shape[0] / self.factor
        ny = plot.image._A.shape[1] / self.factor
        pixX = _MPL.Pixelize(plot.data['px'],
                             plot.data['py'],
                             plot.data['pdx'],
                             plot.data['pdy'],
                             plot.data[self.field_x] - self.bv_x,
                             int(nx), int(ny),
                           (x0, x1, y0, y1),).transpose()
        pixY = _MPL.Pixelize(plot.data['px'],
                             plot.data['py'],
                             plot.data['pdx'],
                             plot.data['pdy'],
                             plot.data[self.field_y] - self.bv_y,
                             int(nx), int(ny),
                           (x0, x1, y0, y1),).transpose()
        X,Y = np.meshgrid(np.linspace(xx0,xx1,nx,endpoint=True),
                          np.linspace(yy0,yy1,ny,endpoint=True))
        if self.normalize:
            nn = np.sqrt(pixX**2 + pixY**2)
            pixX /= nn
            pixY /= nn
        plot._axes.quiver(X,Y, pixX, pixY, scale=self.scale, scale_units=self.scale_units)
        plot._axes.set_xlim(xx0,xx1)
        plot._axes.set_ylim(yy0,yy1)
        plot._axes.hold(False)

class ContourCallback(PlotCallback):
    _type_name = "contour"
    def __init__(self, field, ncont=5, factor=4, clim=None,
                 plot_args = None, label = False, label_args = None):
        """
        annotate_contour(self, field, ncont=5, factor=4, take_log=False, clim=None,
                         plot_args = None):

        Add contours in *field* to the plot.  *ncont* governs the number of
        contours generated, *factor* governs the number of points used in the
        interpolation, *take_log* governs how it is contoured and *clim* gives
        the (upper, lower) limits for contouring.
        """
        PlotCallback.__init__(self)
        self.ncont = ncont
        self.field = field
        self.factor = factor
        from yt.utilities.delaunay.triangulate import Triangulation as triang
        self.triang = triang
        self.clim = clim
        if plot_args is None: plot_args = {'colors':'k'}
        self.plot_args = plot_args
        self.label = label
        if label_args is None:
            label_args = {}
        self.label_args = label_args

    def __call__(self, plot):
        x0, x1 = plot.xlim
        y0, y1 = plot.ylim
        
        xx0, xx1 = plot._axes.get_xlim()
        yy0, yy1 = plot._axes.get_ylim()
        
        plot._axes.hold(True)
        
        numPoints_x = plot.image._A.shape[0]
        numPoints_y = plot.image._A.shape[1]
        
        # Multiply by dx and dy to go from data->plot
        dx = (xx1 - xx0) / (x1-x0)
        dy = (yy1 - yy0) / (y1-y0)

        #dcollins Jan 11 2009.  Improved to allow for periodic shifts in the plot.
        #Now makes a copy of the position fields "px" and "py" and adds the
        #appropriate shift to the coppied field.  

        #set the cumulative arrays for the periodic shifting.
        AllX = np.zeros(plot.data["px"].size, dtype='bool')
        AllY = np.zeros(plot.data["py"].size, dtype='bool')
        XShifted = plot.data["px"].copy()
        YShifted = plot.data["py"].copy()
        dom_x, dom_y = plot._period
        for shift in np.mgrid[-1:1:3j]:
            xlim = ((plot.data["px"] + shift*dom_x >= x0)
                 &  (plot.data["px"] + shift*dom_x <= x1))
            ylim = ((plot.data["py"] + shift*dom_y >= y0)
                 &  (plot.data["py"] + shift*dom_y <= y1))
            XShifted[xlim] += shift * dom_x
            YShifted[ylim] += shift * dom_y
            AllX |= xlim
            AllY |= ylim
        # At this point XShifted and YShifted are the shifted arrays of
        # position data in data coordinates
        wI = (AllX & AllY)

        # We want xi, yi in plot coordinates
        xi, yi = np.mgrid[xx0:xx1:numPoints_x/(self.factor*1j),\
                          yy0:yy1:numPoints_y/(self.factor*1j)]

        # This converts XShifted and YShifted into plot coordinates
        x = (XShifted[wI]-x0)*dx + xx0
        y = (YShifted[wI]-y0)*dy + yy0
        z = plot.data[self.field][wI]
        if plot.pf.field_info[self.field].take_log: z=np.log10(z)

        # Both the input and output from the triangulator are in plot
        # coordinates
        zi = self.triang(x,y).nn_interpolator(z)(xi,yi)
        
        if plot.pf.field_info[self.field].take_log and self.clim is not None: 
            self.clim = (np.log10(self.clim[0]), np.log10(self.clim[1]))
        
        if self.clim is not None: 
            self.ncont = np.linspace(self.clim[0], self.clim[1], self.ncont)
        
        cset = plot._axes.contour(xi,yi,zi,self.ncont, **self.plot_args)
        plot._axes.set_xlim(xx0,xx1)
        plot._axes.set_ylim(yy0,yy1)
        plot._axes.hold(False)
        
        if self.label:
            plot._axes.clabel(cset, **self.label_args)
        

class GridBoundaryCallback(PlotCallback):
    _type_name = "grids"
    def __init__(self, alpha=1.0, min_pix=1, annotate=False, periodic=True):
        """ 
        annotate_grids(alpha=1.0, min_pix=1, annotate=False, periodic=True)

        Adds grid boundaries to a plot, optionally with *alpha*-blending.
        Cuttoff for display is at *min_pix* wide.
        *annotate* puts the grid id in the corner of the grid.  (Not so great in projections...)
        """
        PlotCallback.__init__(self)
        self.alpha = alpha
        self.min_pix = min_pix
        self.annotate = annotate # put grid numbers in the corner.
        self.periodic = periodic

    def __call__(self, plot):
        x0, x1 = plot.xlim
        y0, y1 = plot.ylim
        width, height = plot.image._A.shape
        xx0, xx1 = plot._axes.get_xlim()
        yy0, yy1 = plot._axes.get_ylim()
        xi = x_dict[plot.data.axis]
        yi = y_dict[plot.data.axis]
        dx = width / (x1-x0)
        dy = height / (y1-y0)
        px_index = x_dict[plot.data.axis]
        py_index = y_dict[plot.data.axis]
        dom = plot.data.pf.domain_right_edge - plot.data.pf.domain_left_edge
        if self.periodic:
            pxs, pys = np.mgrid[-1:1:3j,-1:1:3j]
        else:
            pxs, pys = np.mgrid[0:0:1j,0:0:1j]
        GLE = plot.data.pf.h.grid_left_edge
        GRE = plot.data.pf.h.grid_right_edge
        for px_off, py_off in zip(pxs.ravel(), pys.ravel()):
            pxo = px_off * dom[px_index]
            pyo = py_off * dom[py_index]
            left_edge_px = (GLE[:,px_index]+pxo-x0)*dx
            left_edge_py = (GLE[:,py_index]+pyo-y0)*dy
            right_edge_px = (GRE[:,px_index]+pxo-x0)*dx
            right_edge_py = (GRE[:,py_index]+pyo-y0)*dy
            verts = np.array(
                [(left_edge_px, left_edge_px, right_edge_px, right_edge_px),
                 (left_edge_py, right_edge_py, right_edge_py, left_edge_py)])
            visible =  (right_edge_px - left_edge_px > self.min_pix) & \
                       (right_edge_px - left_edge_px > self.min_pix)
            verts=verts.transpose()[visible,:,:]
            if verts.size == 0: continue
            edgecolors = (0.0,0.0,0.0,self.alpha)
            verts[:,:,0]= (xx1-xx0)*(verts[:,:,0]/width) + xx0
            verts[:,:,1]= (yy1-yy0)*(verts[:,:,1]/height) + yy0
            grid_collection = matplotlib.collections.PolyCollection(
                verts, facecolors="none",
                edgecolors=edgecolors)
            plot._axes.hold(True)
            plot._axes.add_collection(grid_collection)
            if self.annotate:
                ids = [g.id for g in plot.data._grids]
                for n in range(len(left_edge_px)):
                    plot._axes.text(left_edge_px[n]+2,left_edge_py[n]+2,ids[n])
            plot._axes.hold(False)

class StreamlineCallback(PlotCallback):
    _type_name = "streamlines"
    def __init__(self, field_x, field_y, factor = 16,
                 density = 1, plot_args=None):
        """
        annotate_streamlines(field_x, field_y, factor = 16,
                             density = 1, plot_args=None):

        Add streamlines to any plot, using the *field_x* and *field_y*
        from the associated data, skipping every *factor* datapoints like
        'quiver'. *density* is the index of the amount of the streamlines.
        """
        PlotCallback.__init__(self)
        self.field_x = field_x
        self.field_y = field_y
        self.bv_x = self.bv_y = 0
        self.factor = factor
        self.dens = density
        if plot_args is None: plot_args = {}
        self.plot_args = plot_args
        
    def __call__(self, plot):
        x0, x1 = plot.xlim
        y0, y1 = plot.ylim
        xx0, xx1 = plot._axes.get_xlim()
        yy0, yy1 = plot._axes.get_ylim()
        plot._axes.hold(True)
        nx = plot.image._A.shape[0] / self.factor
        ny = plot.image._A.shape[1] / self.factor
        pixX = _MPL.Pixelize(plot.data['px'],
                             plot.data['py'],
                             plot.data['pdx'],
                             plot.data['pdy'],
                             plot.data[self.field_x] - self.bv_x,
                             int(nx), int(ny),
                           (x0, x1, y0, y1),).transpose()
        pixY = _MPL.Pixelize(plot.data['px'],
                             plot.data['py'],
                             plot.data['pdx'],
                             plot.data['pdy'],
                             plot.data[self.field_y] - self.bv_y,
                             int(nx), int(ny),
                           (x0, x1, y0, y1),).transpose()
        X,Y = (np.linspace(xx0,xx1,nx,endpoint=True),
                          np.linspace(yy0,yy1,ny,endpoint=True))
        plot._axes.streamplot(X,Y, pixX, pixY, density = self.dens,
                              **self.plot_args)
        plot._axes.set_xlim(xx0,xx1)
        plot._axes.set_ylim(yy0,yy1)
        plot._axes.hold(False)

class LabelCallback(PlotCallback):
    _type_name = "axis_label"
    def __init__(self, label):
        """
        This adds a label to the plot.
        """
        PlotCallback.__init__(self)
        self.label = label

    def __call__(self, plot):
        plot._figure.subplots_adjust(hspace=0, wspace=0,
                                     bottom=0.1, top=0.9,
                                     left=0.0, right=1.0)
        plot._axes.set_xlabel(self.label)
        plot._axes.set_ylabel(self.label)

def get_smallest_appropriate_unit(v, pf):
    max_nu = 1e30
    good_u = None
    for unit in ['mpc','kpc','pc','au','rsun','cm']:
        vv = v*pf[unit]
        if vv < max_nu and vv > 1.0:
            good_u = unit
            max_nu = v*pf[unit]
    if good_u is None : good_u = 'cm'
    return good_u

class UnitBoundaryCallback(PlotCallback):
    _type_name = "units"
    def __init__(self, unit = "au", factor=4, text_annotate=True, text_which=-2):
        """
        Add on a plot indicating where *factor*s of *unit* are shown.
        Optionally *text_annotate* on the *text_which*-indexed box on display.
        """
        PlotCallback.__init__(self)
        self.unit = unit
        self.factor = factor
        self.text_annotate = text_annotate
        self.text_which = -2

    def __call__(self, plot):
        x0, x1 = plot.xlim
        y0, y1 = plot.ylim
        l, b, width, height = mpl_get_bounds(plot._axes.bbox)
        xi = x_dict[plot.data.axis]
        yi = y_dict[plot.data.axis]
        dx = plot.image._A.shape[0] / (x1-x0)
        dy = plot.image._A.shape[1] / (y1-y0)
        center = plot.data.center
        min_dx = plot.data['pdx'].min()
        max_dx = plot.data['pdx'].max()
        w_min_x = 250.0 * min_dx
        w_max_x = 1.0 / self.factor
        min_exp_x = np.ceil(np.log10(w_min_x*plot.data.pf[self.unit])
                           /np.log10(self.factor))
        max_exp_x = np.floor(np.log10(w_max_x*plot.data.pf[self.unit])
                            /np.log10(self.factor))
        n_x = max_exp_x - min_exp_x + 1
        widths = np.logspace(min_exp_x, max_exp_x, num = n_x, base=self.factor)
        widths /= plot.data.pf[self.unit]
        left_edge_px = (center[xi] - widths/2.0 - x0)*dx
        left_edge_py = (center[yi] - widths/2.0 - y0)*dy
        right_edge_px = (center[xi] + widths/2.0 - x0)*dx
        right_edge_py = (center[yi] + widths/2.0 - y0)*dy
        verts = np.array(
                [(left_edge_px, left_edge_px, right_edge_px, right_edge_px),
                 (left_edge_py, right_edge_py, right_edge_py, left_edge_py)])
        visible =  ( right_edge_px - left_edge_px > 25 ) & \
                   ( right_edge_px - left_edge_px > 25 ) & \
                   ( (right_edge_px < width) & (left_edge_px > 0) ) & \
                   ( (right_edge_py < height) & (left_edge_py > 0) )
        verts=verts.transpose()[visible,:,:]
        grid_collection = matplotlib.collections.PolyCollection(
                verts, facecolors="none",
                       edgecolors=(0.0,0.0,0.0,1.0),
                       linewidths=2.5)
        plot._axes.hold(True)
        plot._axes.add_collection(grid_collection)
        if self.text_annotate:
            ti = max(self.text_which, -1*len(widths[visible]))
            if ti < len(widths[visible]): 
                w = widths[visible][ti]
                good_u = get_smallest_appropriate_unit(w, plot.data.pf)
                w *= plot.data.pf[good_u]
                plot._axes.annotate("%0.3e %s" % (w,good_u), verts[ti,1,:]+5)
        plot._axes.hold(False)

class LinePlotCallback(PlotCallback):
    _type_name = "line"
    def __init__(self, x, y, plot_args = None):
        """
        annotate_line(x, y, plot_args = None)

        Over plot *x* and *y* with *plot_args* fed into the plot.
        """
        PlotCallback.__init__(self)
        self.x = x
        self.y = y
        if plot_args is None: plot_args = {}
        self.plot_args = plot_args

    def __call__(self, plot):
        plot._axes.hold(True)
        plot._axes.plot(self.x, self.y, **self.plot_args)
        plot._axes.hold(False)

class ImageLineCallback(LinePlotCallback):
    _type_name = "image_line"

    def __init__(self, p1, p2, data_coords=False, plot_args = None):
        """
        annotate_image_line(p1, p2, data_coords=False, plot_args = None)

        Plot from *p1* to *p2* (image plane coordinates)
        with *plot_args* fed into the plot.
        """
        PlotCallback.__init__(self)
        self.p1 = p1
        self.p2 = p2
        if plot_args is None: plot_args = {}
        self.plot_args = plot_args
        self._ids = []
        self.data_coords = data_coords

    def __call__(self, plot):
        # We manually clear out any previous calls to this callback:
        plot._axes.lines = [l for l in plot._axes.lines if id(l) not in self._ids]
        kwargs = self.plot_args.copy()
        if self.data_coords and len(plot.image._A.shape) == 2:
            p1 = self.convert_to_plot(plot, self.p1)
            p2 = self.convert_to_plot(plot, self.p2)
        else:
            p1, p2 = self.p1, self.p2
            if not self.data_coords:
                kwargs["transform"] = plot._axes.transAxes

        px, py = (p1[0], p2[0]), (p1[1], p2[1])

        # Save state
        xx0, xx1 = plot._axes.get_xlim()
        yy0, yy1 = plot._axes.get_ylim()
        plot._axes.hold(True)
        ii = plot._axes.plot(px, py, **kwargs)
        self._ids.append(id(ii[0]))
        # Reset state
        plot._axes.set_xlim(xx0,xx1)
        plot._axes.set_ylim(yy0,yy1)
        plot._axes.hold(False)

class CuttingQuiverCallback(PlotCallback):
    _type_name = "cquiver"
    def __init__(self, field_x, field_y, factor):
        """
        annotate_cquiver(field_x, field_y, factor)

        Get a quiver plot on top of a cutting plane, using *field_x* and
        *field_y*, skipping every *factor* datapoint in the discretization.
        """
        PlotCallback.__init__(self)
        self.field_x = field_x
        self.field_y = field_y
        self.factor = factor

    def __call__(self, plot):
        x0, x1 = plot.xlim
        y0, y1 = plot.ylim
        xx0, xx1 = plot._axes.get_xlim()
        yy0, yy1 = plot._axes.get_ylim()
        plot._axes.hold(True)
        nx = plot.image._A.shape[0] / self.factor
        ny = plot.image._A.shape[1] / self.factor
        indices = np.argsort(plot.data['dx'])[::-1]
        pixX = _MPL.CPixelize( plot.data['x'], plot.data['y'], plot.data['z'],
                               plot.data['px'], plot.data['py'],
                               plot.data['pdx'], plot.data['pdy'], plot.data['pdz'],
                               plot.data.center, plot.data._inv_mat, indices,
                               plot.data[self.field_x],
                               int(nx), int(ny),
                               (x0, x1, y0, y1),).transpose()
        pixY = _MPL.CPixelize( plot.data['x'], plot.data['y'], plot.data['z'],
                               plot.data['px'], plot.data['py'],
                               plot.data['pdx'], plot.data['pdy'], plot.data['pdz'],
                               plot.data.center, plot.data._inv_mat, indices,
                               plot.data[self.field_y],
                               int(nx), int(ny),
                               (x0, x1, y0, y1),).transpose()
        X,Y = np.meshgrid(np.linspace(xx0,xx1,nx,endpoint=True),
                          np.linspace(yy0,yy1,ny,endpoint=True))
        plot._axes.quiver(X,Y, pixX, pixY)
        plot._axes.set_xlim(xx0,xx1)
        plot._axes.set_ylim(yy0,yy1)
        plot._axes.hold(False)

class ClumpContourCallback(PlotCallback):
    _type_name = "clumps"
    def __init__(self, clumps, plot_args = None):
        """
        annotate_clumps(clumps, plot_args = None)

        Take a list of *clumps* and plot them as a set of contours.
        """
        self.clumps = clumps
        if plot_args is None: plot_args = {}
        self.plot_args = plot_args

    def __call__(self, plot):
        x0, x1 = plot.xlim
        y0, y1 = plot.ylim
        xx0, xx1 = plot._axes.get_xlim()
        yy0, yy1 = plot._axes.get_ylim()
        plot._axes.hold(True)

        px_index = x_dict[plot.data.axis]
        py_index = y_dict[plot.data.axis]

        xf = axis_names[px_index]
        yf = axis_names[py_index]

        DomainRight = plot.data.pf.domain_right_edge
        DomainLeft = plot.data.pf.domain_left_edge
        DomainWidth = DomainRight - DomainLeft
        
        nx, ny = plot.image._A.shape
        buff = np.zeros((nx,ny),dtype='float64')
        for i,clump in enumerate(reversed(self.clumps)):
            mylog.debug("Pixelizing contour %s", i)


            xf_copy = clump[xf].copy()
            yf_copy = clump[yf].copy()
            
            temp = _MPL.Pixelize(xf_copy, yf_copy, 
                                 clump['dx']/2.0,
                                 clump['dy']/2.0,
                                 clump['dx']*0.0+i+1, # inits inside Pixelize
                                 int(nx), int(ny),
                             (x0, x1, y0, y1), 0).transpose()
            buff = np.maximum(temp, buff)
        self.rv = plot._axes.contour(buff, len(self.clumps)+1,
                                     **self.plot_args)
        plot._axes.hold(False)

class ArrowCallback(PlotCallback):
    _type_name = "arrow"
    def __init__(self, pos, code_size, plot_args = None):
        """
        annotate_arrow(pos, code_size, plot_args = None)

        This adds an arrow pointing at *pos* with size *code_size* in code
        units.  *plot_args* is a dict fed to matplotlib with arrow properties.
        """
        self.pos = pos
        if not iterable(code_size):
            code_size = (code_size, code_size)
        self.code_size = code_size
        if plot_args is None: plot_args = {}
        self.plot_args = plot_args

    def __call__(self, plot):
        if len(self.pos) == 3:
            pos = (self.pos[x_dict[plot.data.axis]],
                   self.pos[y_dict[plot.data.axis]])
        else: pos = self.pos
        from matplotlib.patches import Arrow
        # Now convert the pixels to code information
        x, y = self.convert_to_plot(plot, pos)
        dx, dy = self.convert_to_plot(plot, self.code_size, False)
        arrow = Arrow(x, y, dx, dy, **self.plot_args)
        plot._axes.add_patch(arrow)

class PointAnnotateCallback(PlotCallback):
    _type_name = "point"
    def __init__(self, pos, text, text_args = None):
        """
        annotate_point(pos, text, text_args = None)

        This adds *text* at position *pos*, where *pos* is in code-space.
        *text_args* is a dict fed to the text placement code.
        """
        self.pos = pos
        self.text = text
        if text_args is None: text_args = {}
        self.text_args = text_args

    def __call__(self, plot):
        if len(self.pos) == 3:
            pos = (self.pos[x_dict[plot.data.axis]],
                   self.pos[y_dict[plot.data.axis]])
        else: pos = self.pos
        width,height = plot.image._A.shape
        x,y = self.convert_to_plot(plot, pos)
        
        plot._axes.text(x, y, self.text, **self.text_args)

class MarkerAnnotateCallback(PlotCallback):
    _type_name = "marker"
    def __init__(self, pos, marker='x', plot_args=None):
        """
        annotate_marker(pos, marker='x', plot_args=None)

        Adds text *marker* at *pos* in code-arguments.  *plot_args* is a dict
        that will be forwarded to the plot command.
        """
        self.pos = pos
        self.marker = marker
        if plot_args is None: plot_args = {}
        self.plot_args = plot_args

    def __call__(self, plot):
        xx0, xx1 = plot._axes.get_xlim()
        yy0, yy1 = plot._axes.get_ylim()
        if np.array(self.pos).shape == (3,):
            pos = (self.pos[x_dict[plot.data.axis]],
                   self.pos[y_dict[plot.data.axis]])
        elif np.array(self.pos).shape == (2,):
            pos = self.pos
        x,y = self.convert_to_plot(plot, pos)
        plot._axes.hold(True)
        plot._axes.scatter(x,y, marker = self.marker, **self.plot_args)
        plot._axes.set_xlim(xx0,xx1)
        plot._axes.set_ylim(yy0,yy1)
        plot._axes.hold(False)

class SphereCallback(PlotCallback):
    _type_name = "sphere"
    def __init__(self, center, radius, circle_args = None,
                 text = None, text_args = None):
        """
        annotate_sphere(center, radius, circle_args = None,
                        text = None, text_args = None)
        
        A sphere centered at *center* in code units with radius *radius* in
        code units will be created, with optional *circle_args*, *text*, and
        *text_args*.
        """
        self.center = center
        self.radius = radius
        if circle_args is None: circle_args = {}
        if 'fill' not in circle_args: circle_args['fill'] = False
        self.circle_args = circle_args
        self.text = text
        self.text_args = text_args
        if self.text_args is None: self.text_args = {}

    def __call__(self, plot):
        from matplotlib.patches import Circle
        
        radius = self.radius * self.pixel_scale(plot)[0]

        (xi, yi) = (x_dict[plot.data.axis], y_dict[plot.data.axis])

        (center_x,center_y) = self.convert_to_plot(plot,(self.center[xi], self.center[yi]))
        
        cir = Circle((center_x, center_y), radius, **self.circle_args)
        plot._axes.add_patch(cir)
        if self.text is not None:
            plot._axes.text(center_x, center_y, self.text,
                            **self.text_args)

class HopCircleCallback(PlotCallback):
    _type_name = "hop_circles"
    def __init__(self, hop_output, max_number=None,
                 annotate=False, min_size=20, max_size=10000000,
                 font_size=8, print_halo_size=False,
                 print_halo_mass=False, width=None):
        """
        annotate_hop_circles(hop_output, max_number=None,
                             annotate=False, min_size=20, max_size=10000000,
                             font_size=8, print_halo_size=False,
                             print_halo_mass=False, width=None)

        Accepts a :class:`yt.HopList` *hop_output* and plots up to
        *max_number* (None for unlimited) halos as circles.
        """
        self.hop_output = hop_output
        self.max_number = max_number
        self.annotate = annotate
        self.min_size = min_size
        self.max_size = max_size
        self.font_size = font_size
        self.print_halo_size = print_halo_size
        self.print_halo_mass = print_halo_mass
        self.width = width

    def __call__(self, plot):
        from matplotlib.patches import Circle
        for halo in self.hop_output[:self.max_number]:
            size = halo.get_size()
            if size < self.min_size or size > self.max_size: continue
            # This could use halo.maximum_radius() instead of width
            if self.width is not None and \
                np.abs(halo.center_of_mass() - 
                       plot.data.center)[plot.data.axis] > \
                   self.width:
                continue
            
            radius = halo.maximum_radius() * self.pixel_scale(plot)[0]
            center = halo.center_of_mass()
            
            (xi, yi) = (x_dict[plot.data.axis], y_dict[plot.data.axis])

            (center_x,center_y) = self.convert_to_plot(plot,(center[xi], center[yi]))
            cir = Circle((center_x, center_y), radius, fill=False)
            plot._axes.add_patch(cir)
            if self.annotate:
                if self.print_halo_size:
                    plot._axes.text(center_x, center_y, "%s" % size,
                    fontsize=self.font_size)
                elif self.print_halo_mass:
                    plot._axes.text(center_x, center_y, "%s" % halo.total_mass(),
                    fontsize=self.font_size)
                else:
                    plot._axes.text(center_x, center_y, "%s" % halo.id,
                    fontsize=self.font_size)

class HopParticleCallback(PlotCallback):
    _type_name = "hop_particles"
    def __init__(self, hop_output, max_number=None, p_size=1.0,
                min_size=20, alpha=0.2):
        """
        annotate_hop_particles(hop_output, max_number, p_size=1.0,
                               min_size=20, alpha=0.2):

        Adds particle positions for the members of each halo as identified
        by HOP. Along *axis* up to *max_number* groups in *hop_output* that are
        larger than *min_size* are plotted with *p_size* pixels per particle; 
        *alpha* determines the opacity of each particle.
        """
        self.hop_output = hop_output
        self.p_size = p_size
        self.max_number = max_number
        self.min_size = min_size
        self.alpha = alpha
    
    def __call__(self,plot):
        (dx,dy) = self.pixel_scale(plot)

        (xi, yi) = (x_names[plot.data.axis], y_names[plot.data.axis])

        # now we loop over the haloes
        for halo in self.hop_output[:self.max_number]:
            size = halo.get_size()

            if size < self.min_size: continue

            (px,py) = self.convert_to_plot(plot,(halo["particle_position_%s" % xi],
                                                 halo["particle_position_%s" % yi]))
            
            # Need to get the plot limits and set the hold state before scatter
            # and then restore the limits and turn off the hold state afterwards
            # because scatter will automatically adjust the plot window which we
            # do not want
            
            xlim = plot._axes.get_xlim()
            ylim = plot._axes.get_ylim()
            plot._axes.hold(True)

            plot._axes.scatter(px, py, edgecolors="None",
                s=self.p_size, c='black', alpha=self.alpha)
            
            plot._axes.set_xlim(xlim)
            plot._axes.set_ylim(ylim)
            plot._axes.hold(False)


class CoordAxesCallback(PlotCallback):
    _type_name = "coord_axes"
    def __init__(self,unit=None,coords=False):
        """
        Creates x and y axes for a VMPlot. In the future, it will
        attempt to guess the proper units to use.
        """
        PlotCallback.__init__(self)
        self.unit = unit
        self.coords = coords

    def __call__(self,plot):
        # 1. find out what the domain is
        # 2. pick a unit for it
        # 3. run self._axes.set_xlabel & self._axes.set_ylabel to actually lay things down.
        # 4. adjust extent information to make sure labels are visable.

        # put the plot into data coordinates
        nx,ny = plot.image._A.shape
        dx = (plot.xlim[1] - plot.xlim[0])/nx
        dy = (plot.ylim[1] - plot.ylim[0])/ny

        unit_conversion = plot.pf[plot.im["Unit"]]
        aspect = (plot.xlim[1]-plot.xlim[0])/(plot.ylim[1]-plot.ylim[0])

        print "aspect ratio = ", aspect

        # if coords is False, label axes relative to the center of the
        # display. if coords is True, label axes with the absolute
        # coordinates of the region.
        xcenter = 0.
        ycenter = 0.
        if not self.coords:
            center = plot.data.center
            if plot.data.axis == 0:
                xcenter = center[1]
                ycenter = center[2]
            elif plot.data.axis == 1:
                xcenter = center[0]
                ycenter = center[2]
            else:
                xcenter = center[0]
                ycenter = center[1]


            xformat_function = lambda a,b: '%7.1e' %((a*dx + plot.xlim[0] - xcenter)*unit_conversion)
            yformat_function = lambda a,b: '%7.1e' %((a*dy + plot.ylim[0] - ycenter)*unit_conversion)
        else:
            xformat_function = lambda a,b: '%7.1e' %((a*dx + plot.xlim[0])*unit_conversion)
            yformat_function = lambda a,b: '%7.1e' %((a*dy + plot.ylim[0])*unit_conversion)
            
        xticker = matplotlib.ticker.FuncFormatter(xformat_function)
        yticker = matplotlib.ticker.FuncFormatter(yformat_function)
        plot._axes.xaxis.set_major_formatter(xticker)
        plot._axes.yaxis.set_major_formatter(yticker)
        
        xlabel = '%s (%s)' % (axis_labels[plot.data.axis][0],plot.im["Unit"])
        ylabel = '%s (%s)' % (axis_labels[plot.data.axis][1],plot.im["Unit"])
        xticksize = nx/4.
        yticksize = ny/4.
        plot._axes.xaxis.set_major_locator(matplotlib.ticker.FixedLocator([i*xticksize for i in range(0,5)]))
        plot._axes.yaxis.set_major_locator(matplotlib.ticker.FixedLocator([i*yticksize for i in range(0,5)]))
        
        plot._axes.set_xlabel(xlabel,visible=True)
        plot._axes.set_ylabel(ylabel,visible=True)
        plot._figure.subplots_adjust(left=0.1,right=0.8)

class TextLabelCallback(PlotCallback):
    _type_name = "text"
    def __init__(self, pos, text, data_coords=False, text_args = None):
        """
        annotate_text(pos, text, data_coords=False, text_args = None)

        Accepts a position in (0..1, 0..1) of the image, some text and
        optionally some text arguments. If data_coords is True,
        position will be in code units instead of image coordinates.
        """
        self.pos = pos
        self.text = text
        self.data_coords = data_coords
        if text_args is None: text_args = {}
        self.text_args = text_args

    def __call__(self, plot):
        kwargs = self.text_args.copy()
        if self.data_coords and len(plot.image._A.shape) == 2:
            if len(self.pos) == 3:
                pos = (self.pos[x_dict[plot.data.axis]],
                       self.pos[y_dict[plot.data.axis]])
            else: pos = self.pos
            x,y = self.convert_to_plot(plot, pos)
        else:
            x, y = self.pos
            if not self.data_coords:
                kwargs["transform"] = plot._axes.transAxes
        plot._axes.text(x, y, self.text, **kwargs)

class ParticleCallback(PlotCallback):
    _type_name = "particles"
    region = None
    _descriptor = None
    def __init__(self, width, p_size=1.0, col='k', marker='o', stride=1.0,
                 ptype=None, stars_only=False, dm_only=False,
                 minimum_mass=None, alpha=1.0):
        """
        annotate_particles(width, p_size=1.0, col='k', marker='o', stride=1.0,
                           ptype=None, stars_only=False, dm_only=False,
                           minimum_mass=None, alpha=1.0)

        Adds particle positions, based on a thick slab along *axis* with a
        *width* along the line of sight.  *p_size* controls the number of
        pixels per particle, and *col* governs the color.  *ptype* will
        restrict plotted particles to only those that are of a given type.
        *minimum_mass* will require that the particles be of a given mass,
        calculated via ParticleMassMsun, to be plotted. *alpha* determines
        each particle's opacity.
        """
        PlotCallback.__init__(self)
        self.width = width
        self.p_size = p_size
        self.color = col
        self.marker = marker
        self.stride = stride
        self.ptype = ptype
        self.stars_only = stars_only
        self.dm_only = dm_only
        self.minimum_mass = minimum_mass
        self.alpha = alpha

    def __call__(self, plot):
        data = plot.data
        # we construct a recantangular prism
        x0, x1 = plot.xlim
        y0, y1 = plot.ylim
        xx0, xx1 = plot._axes.get_xlim()
        yy0, yy1 = plot._axes.get_ylim()
        reg = self._get_region((x0,x1), (y0,y1), plot.data.axis, data)
        field_x = "particle_position_%s" % axis_names[x_dict[data.axis]]
        field_y = "particle_position_%s" % axis_names[y_dict[data.axis]]
        gg = ( ( reg[field_x] >= x0 ) & ( reg[field_x] <= x1 )
           &   ( reg[field_y] >= y0 ) & ( reg[field_y] <= y1 ) )
        if self.ptype is not None:
            gg &= (reg["particle_type"] == self.ptype)
            if gg.sum() == 0: return
        if self.stars_only:
            gg &= (reg["creation_time"] > 0.0)
            if gg.sum() == 0: return
        if self.dm_only:
            gg &= (reg["creation_time"] <= 0.0)
            if gg.sum() == 0: return
        if self.minimum_mass is not None:
            gg &= (reg["ParticleMassMsun"] >= self.minimum_mass)
            if gg.sum() == 0: return
        plot._axes.hold(True)
        px, py = self.convert_to_plot(plot,
                    [reg[field_x][gg][::self.stride],
                     reg[field_y][gg][::self.stride]])
        plot._axes.scatter(px, py, edgecolors='None', marker=self.marker,
                           s=self.p_size, c=self.color,alpha=self.alpha)
        plot._axes.set_xlim(xx0,xx1)
        plot._axes.set_ylim(yy0,yy1)
        plot._axes.hold(False)


    def _get_region(self, xlim, ylim, axis, data):
        LE, RE = [None]*3, [None]*3
        xax = x_dict[axis]
        yax = y_dict[axis]
        zax = axis
        LE[xax], RE[xax] = xlim
        LE[yax], RE[yax] = ylim
        LE[zax] = data.center[zax] - self.width*0.5
        RE[zax] = data.center[zax] + self.width*0.5
        if self.region is not None \
            and np.all(self.region.left_edge <= LE) \
            and np.all(self.region.right_edge >= RE):
            return self.region
        self.region = data.pf.h.periodic_region(
            data.center, LE, RE)
        return self.region

class TitleCallback(PlotCallback):
    _type_name = "title"
    def __init__(self, title):
        """
        annotate_title(title)

        Accepts a *title* and adds it to the plot
        """
        PlotCallback.__init__(self)
        self.title = title

    def __call__(self,plot):
        plot._axes.set_title(self.title)

class FlashRayDataCallback(PlotCallback):
    _type_name = "flash_ray_data"
    def __init__(self, cmap_name='bone', sample=None):
        """ 
        annotate_flash_ray_data(cmap_name='bone', sample=None)

        Adds ray trace data to the plot.  *cmap_name* is the name of the color map 
        ('bone', 'jet', 'hot', etc).  *sample* dictates the amount of down sampling 
        to do to prevent all of the rays from being  plotted.  This may be None 
        (plot all rays, default), an integer (step size), or a slice object.
        """
        self.cmap_name = cmap_name
        self.sample = sample if isinstance(sample, slice) else slice(None, None, sample)

    def __call__(self, plot):
        ray_data = plot.data.pf._handle["RayData"][:]
        idx = ray_data[:,0].argsort(kind="mergesort")
        ray_data = ray_data[idx]

        tags = ray_data[:,0]
        coords = ray_data[:,1:3]
        power = ray_data[:,4]
        power /= power.max()
        cx, cy = self.convert_to_plot(plot, coords.T)
        coords[:,0], coords[:,1] = cx, cy
        splitidx = np.argwhere(0 < (tags[1:] - tags[:-1])) + 1
        coords = np.split(coords, splitidx.flat)[self.sample]
        power = np.split(power, splitidx.flat)[self.sample]
        cmap = matplotlib.cm.get_cmap(self.cmap_name)

        plot._axes.hold(True)
        colors = [cmap(p.max()) for p in power]
        lc = matplotlib.collections.LineCollection(coords, colors=colors)
        plot._axes.add_collection(lc)
        plot._axes.hold(False)


class TimestampCallback(PlotCallback):
    _type_name = "timestamp"
    _time_conv = {
          'as': 1e-18,
          'attosec': 1e-18,
          'attosecond': 1e-18,
          'attoseconds': 1e-18,
          'fs': 1e-15,
          'femtosec': 1e-15,
          'femtosecond': 1e-15,
          'femtoseconds': 1e-15,
          'ps': 1e-12,
          'picosec': 1e-12,
          'picosecond': 1e-12,
          'picoseconds': 1e-12,
          'ns': 1e-9,
          'nanosec': 1e-9,
          'nanosecond':1e-9,
          'nanoseconds' : 1e-9,
          'us': 1e-6,
          'microsec': 1e-6,
          'microsecond': 1e-6,
          'microseconds': 1e-6,
          'ms': 1e-3,
          'millisec': 1e-3,
          'millisecond': 1e-3,
          'milliseconds': 1e-3,
          's': 1.0,
          'sec': 1.0,
          'second':1.0,
          'seconds': 1.0,
          'm': 60.0,
          'min': 60.0,
          'minute': 60.0,
          'minutes': 60.0,
          'h': 3600.0,
          'hour': 3600.0,
          'hours': 3600.0,
          'd': 86400.0,
          'day': 86400.0,
          'days': 86400.0,
          'y': 86400.0*365.25,
          'year': 86400.0*365.25,
          'years': 86400.0*365.25,
          'ev': 1e-9 * 7.6e-8 / 6.03,
          'kev': 1e-12 * 7.6e-8 / 6.03,
          'mev': 1e-15 * 7.6e-8 / 6.03,
          }
    _bbox_dict = {'boxstyle': 'square,pad=0.6', 'fc': 'white', 'ec': 'black', 'alpha': 1.0}

    def __init__(self, x, y, units=None, format="{time:.3G} {units}", normalized = False, 
                 bbox_dict = None, **kwargs):
        """ 
        annotate_timestamp(x, y, units=None, format="{time:.3G} {units}", **kwargs)

        Adds the current time to the plot at point given by *x* and *y*.  If *units* 
        is given ('s', 'ms', 'ns', etc), it will covert the time to this basis.  If 
        *units* is None, it will attempt to figure out the correct value by which to 
        scale.  The *format* keyword is a template string that will be evaluated and 
        displayed on the plot.  If *normalized* is true, *x* and *y* are interpreted 
        as normalized plot coordinates (0,0 is lower-left and 1,1 is upper-right) 
        otherwise *x* and *y* are assumed to be in plot coordinates. The *bbox_dict* 
        is an optional dict of arguments for the bbox that frames the timestamp, see 
        matplotlib's text annotation guide for more details. All other *kwargs* will 
        be passed to the text() method on the plot axes.  See matplotlib's text() 
        functions for more information.
        """
        self.x = x
        self.y = y
        self.format = format
        self.units = units
        self.normalized = normalized
        if bbox_dict is not None:
            self.bbox_dict = bbox_dict
        else:
            self.bbox_dict = self._bbox_dict
        self.kwargs = {'color': 'w'}
        self.kwargs.update(kwargs)

    def __call__(self, plot):
        if self.units is None:
            t = plot.data.pf.current_time * plot.data.pf['Time']
            scale_keys = ['as', 'fs', 'ps', 'ns', 'us', 'ms', 's']
            self.units = 's'
            for k in scale_keys:
                if t < self._time_conv[k]:
                    break
                self.units = k
        t = plot.data.pf.current_time * plot.data.pf['Time'] 
        t /= self._time_conv[self.units.lower()]
        if self.units == 'us':
            self.units = '$\\mu s$'
        s = self.format.format(time=t, units=self.units)
        plot._axes.hold(True)
        if self.normalized:
            plot._axes.text(self.x, self.y, s, horizontalalignment='center',
                            verticalalignment='center', 
                            transform = plot._axes.transAxes, bbox=self.bbox_dict)
        else:
            plot._axes.text(self.x, self.y, s, bbox=self.bbox_dict, **self.kwargs)
        plot._axes.hold(False)


class MaterialBoundaryCallback(ContourCallback):
    _type_name = "material_boundary"
    def __init__(self, field='targ', ncont=1, factor=4, clim=(0.9, 1.0), **kwargs):
        """ 
        annotate_material_boundary(self, field='targ', ncont=1, factor=4, 
                                   clim=(0.9, 1.0), **kwargs):

        Add the limiting contours of *field* to the plot.  Nominally, *field* is 
        the target material but may be any other field present in the hierarchy.
        The number of contours generated is given by *ncount*, *factor* governs 
        the number of points used in the interpolation, and *clim* gives the 
        (upper, lower) limits for contouring.  For this to truly be the boundary
        *clim* should be close to the edge.  For example the default is (0.9, 1.0)
        for 'targ' which is defined on the range [0.0, 1.0].  All other *kwargs* 
        will be passed to the contour() method on the plot axes.  See matplotlib
        for more information.
        """
        plot_args = {'colors': 'w'}
        plot_args.update(kwargs)
        super(MaterialBoundaryCallback, self).__init__(field=field, ncont=ncont, 
                                                       factor=factor, clim=clim,
                                                       plot_args=plot_args)

    def __call__(self, plot):
        super(MaterialBoundaryCallback, self).__call__(plot)