-
Notifications
You must be signed in to change notification settings - Fork 280
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
chrishavlin
wants to merge
2
commits into
yt-project:main
Choose a base branch
from
chrishavlin:fix_depth_handling_off_axis
base: main
Could not load branches
Branch not found: {{ refName }}
Loading
Could not load tags
Nothing to show
Loading
Are you sure you want to change the base?
Some commits from the old base branch may be removed from the timeline,
and old review comments may become outdated.
Open
Changes from 1 commit
Commits
Show all changes
2 commits
Select commit
Hold shift + click to select a range
File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -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) | ||
|
@@ -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: | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe 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): | ||
|
@@ -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 | ||
|
@@ -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.sqrt(np.sum(ds.domain_width.in_units("code_length") ** 2)) * 1.02 | ||
chrishavlin marked this conversation as resolved.
Show resolved
Hide resolved
|
||
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 | ||
|
@@ -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 | ||
|
Oops, something went wrong.
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
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.