Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

BUG: calculate default depth for off axis projections #5081

Open
wants to merge 2 commits into
base: main
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
3 changes: 1 addition & 2 deletions yt/visualization/fixed_resolution.py
Original file line number Diff line number Diff line change
Expand Up @@ -639,7 +639,6 @@ def _generate_image_and_mask(self, item) -> None:
self.bounds[5] - self.bounds[4],
)
)
depth = dd.depth[0] if dd.depth is not None else None
Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

removing this because the depth sanitization before dd creation takes care of this now.

buff = off_axis_projection(
dd.dd,
dd.center,
Expand All @@ -652,7 +651,7 @@ def _generate_image_and_mask(self, item) -> None:
no_ghost=dd.no_ghost,
interpolated=dd.interpolated,
north_vector=dd.north_vector,
depth=depth,
depth=dd.depth,
method=dd.method,
)
if self.data_source.moment == 2:
Expand Down
27 changes: 15 additions & 12 deletions yt/visualization/plot_window.py
Original file line number Diff line number Diff line change
Expand Up @@ -83,7 +83,11 @@ def get_window_parameters(axis, center, width, ds):


def get_oblique_window_parameters(
normal, center, width, ds, depth=None, get3bounds=False
normal,
center,
width,
ds,
depth=None,
):
center, display_center = ds.coordinates.sanitize_center(center, axis=None)
width = ds.coordinates.sanitize_width(normal, width, depth)
Expand All @@ -101,15 +105,7 @@ def get_oblique_window_parameters(

w = tuple(el.in_units("code_length") for el in width)
bounds = tuple(((2 * (i % 2)) - 1) * w[i // 2] / 2 for i in range(len(w) * 2))
if get3bounds and depth is None:
Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

if we initially calculate a depth equal to the diagonal + a bit, we don't need this extra logic here as the off-axis projection calls to this function will always have a depth.

# off-axis projection, depth not specified
# -> set 'large enough' depth using half the box diagonal + margin
d2 = ds.domain_width[0].in_units("code_length") ** 2
d2 += ds.domain_width[1].in_units("code_length") ** 2
d2 += ds.domain_width[2].in_units("code_length") ** 2
diag = np.sqrt(d2)
bounds = bounds + (-0.51 * diag, 0.51 * diag)
return (bounds, center)
return bounds, center


def get_axes_unit(width, ds):
Expand Down Expand Up @@ -2389,7 +2385,8 @@ class OffAxisProjectionPlot(ProjectionPlot, PWViewerMPL):
depth : A tuple or a float
A tuple containing the depth to project through and the string
key of the unit: (width, 'unit'). If set to a float, code units
are assumed
are assumed. In not set, then a depth equal to the diagonal of
the domain width plus a small margin will be used.
weight_field : string
The name of the weighting field. Set to None for no weight.
max_level: int
Expand Down Expand Up @@ -2466,6 +2463,13 @@ def __init__(
"currently supported geometries:"
f" {self._supported_geometries!r}"
)

if depth is None:
# off-axis projection, depth not specified
# -> set 'large enough' depth using half the box diagonal + margin
depth = np.linalg.norm(ds.domain_width.in_units("code_length")) * 1.02
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Why 1.02? Does it need to be an epsilon that large?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Why 1.02?

Good question -- I was just preserving the 2% epsilon in this line

bounds = bounds + (-0.51 * diag, 0.51 * diag)

Does it need to be an epsilon that large?

I suspect it could be smaller, and for periodic datasets in particular a 2% epsilon might not be ideal... @nastasha-w was there reasoning behind the magnitude of margin you used?

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

At least for SPH data, periodicities are applied before rotation, and are never taken to be periodic in the rotated coordinates. Therefore, the large margin shouldn't cause trouble, but a smaller value ought to be fine as well.

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I guess I also figured that the half box diagonal would be a very large margin for most line of sight directions anyway

depth = ds.coordinates.sanitize_depth(depth)[0]

# center_rot normalizes the center to (0,0),
# units match bounds
# for SPH data, we want to input the original center
Expand All @@ -2479,7 +2483,6 @@ def __init__(
width,
ds,
depth=depth,
get3bounds=True,
)
# will probably fail if you try to project an SPH and non-SPH
# field in a single call
Expand Down
Loading