Hi Mark,

This is caused by the fact that OctreeSubsetBlockSlicePosition is not a subclass of YTSelectionContainer, which implements the appropriate API.

I don't remember offhand why it's designed that way. The fix would either be to make OctreeSubsetBlockSlicePosition inherit from YTSelectionContainer (no idea how hard this would be or what it would break to do that) or simply implement the _field_parameter_state context manager for the OctreeSubsetBlockSlicePosition class.

For reference, here's the implementation of the field_parameter_state context manager for the base data object class:

https://github.com/yt-project/yt/blob/master/yt/data_objects/data_containers.py#L1084

If you want to take this on to fix the issue you're running into then that's great. If not then I'd be happy to take a look at this this week. For the latter option, it would be great if you could file an issue for this, along with a script that triggers this problem and makes use of one of the public RAMSES test datasets we have on yt-project.org/data. That will make it easier for me to start on the fix.

-Nathan


On Mon, Feb 12, 2018 at 4:56 PM Mark Richardson <mark.richardson.work@gmail.com> wrote:
Hey, I've previously been able to make off-axis projections of ramses datasets with geometric data containers. However, I've recently tried to do this with a cut_region for dens gas, without success. Here is a script that should work with a ramses cosmo dataset:

import yt
import numpy as np

num = 1
Look = [1.,1.,0]
Up = [0.,1.,1.]


c = np.array([0.5,0.5,0.5])
rad = 0.2
width = 2*rad
extra = 1.05
left  = c - rad*extra
right = c + rad*extra
bb = [left,right]

fn = "../output_{0:05d}/info_{0:05d}.txt"
ds = yt.load(fn.format(num),bbox=bb)

region = ds.region(c,left,right)

d20 = region.cut_region(['obj["density"] > 1e-20'])

p1 = yt.OffAxisProjectionPlot(ds,Look,"density", weight_field='density',center=c,width=width,data_source=region,north_vector=Up)
p1.save("OAProjDens_Region_{0:03d}.png".format(nout))

p1 = yt.OffAxisProjectionPlot(ds,Look,"density", weight_field='density',center=c,width=width,data_source=d20,north_vector=Up)
p1.save("OAProjDens_D20_{0:03d}.png".format(nout))


The first OA Projection is generated correctly, but the second one crashes with the message:

yt : [INFO     ] 2018-02-12 22:20:08,915 Making a fixed resolution buffer of (('gas', 'density')) 800 by 800
yt : [INFO     ] 2018-02-12 22:22:05,737 Saving plot OAProjDens_Region_001.png
yt : [INFO     ] 2018-02-12 22:22:06,641 xlim = -0.2 0.2
yt : [INFO     ] 2018-02-12 22:22:06,641 ylim = -0.2 0.2
yt : [INFO     ] 2018-02-12 22:22:06,641 zlim = -0.500000 0.500000
yt : [INFO     ] 2018-02-12 22:23:58,431 Making a fixed resolution buffer of (('gas', 'density')) 800 by 800
Traceback (most recent call last):
  File "MakeOffAxisProjClouds.py", line 176, in <module>
    p1 = yt.OffAxisProjectionPlot(ds,Look,"density", weight_field='density',center=c,width=width,data_source=d20,north_vector=Up)
  File "/mnt/blackwhale/home/mrichard/soft/yt/yt_cphyc/yt/visualization/plot_window.py", line 1727, in __init__
    fontsize=fontsize)
  File "/mnt/blackwhale/home/mrichard/soft/yt/yt_cphyc/yt/visualization/plot_window.py", line 660, in __init__
    PlotWindow.__init__(self, *args, **kwargs)
  File "/mnt/blackwhale/home/mrichard/soft/yt/yt_cphyc/yt/visualization/plot_window.py", line 219, in __init__
    self._setup_plots()
  File "/mnt/blackwhale/home/mrichard/soft/yt/yt_cphyc/yt/visualization/plot_window.py", line 748, in _setup_plots
    self._recreate_frb()
  File "/mnt/blackwhale/home/mrichard/soft/yt/yt_cphyc/yt/visualization/plot_window.py", line 1733, in _recreate_frb
    super(OffAxisProjectionPlot, self)._recreate_frb()
  File "/mnt/blackwhale/home/mrichard/soft/yt/yt_cphyc/yt/visualization/plot_window.py", line 272, in _recreate_frb
    self._frb._get_data_source_fields()
  File "/mnt/blackwhale/home/mrichard/soft/yt/yt_cphyc/yt/visualization/fixed_resolution.py", line 155, in _get_data_source_fields
    self[f]
  File "/mnt/blackwhale/home/mrichard/soft/yt/yt_cphyc/yt/visualization/fixed_resolution.py", line 543, in __getitem__
    method=dd.method)
  File "/mnt/blackwhale/home/mrichard/soft/yt/yt_cphyc/yt/visualization/volume_rendering/off_axis_projection.py", line 189, in off_axis_projection
    for i, (grid, mask) in enumerate(data_source.blocks):
  File "/mnt/blackwhale/home/mrichard/soft/yt/yt_cphyc/yt/data_objects/selection_data_containers.py", line 836, in blocks
    with obj._field_parameter_state(self.field_parameters):
AttributeError: 'OctreeSubsetBlockSlicePosition' object has no attribute '_field_parameter_state'


I'm trying to find how to have d20 inherit _field_parameter_state from region, but thought someone might be able to point me in the right direction a bit faster.

Thanks!
   -Mark

--

Mark Richardson
MAT Postdoctoral Fellow
Department of Astrophysics
_______________________________________________
yt-dev mailing list -- yt-dev@python.org
To unsubscribe send an email to yt-dev-leave@python.org
_______________________________________________
yt-dev mailing list -- yt-dev@python.org
To unsubscribe send an email to yt-dev-leave@python.org