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

Output field coarsening #580

Merged

Conversation

AlexanderSinn
Copy link
Member

@AlexanderSinn AlexanderSinn commented Jul 30, 2021

This PR adds the option the output a lower resolution field instead of the full resolution one used for calculation. For 3D output, this can greatly reduce the output file size, the amount of GPU memory needed, the time it takes to write the diagnostics. Use:

diagnostic.coarsening = 2 2 2

to set the coarsening ratio for x, y and z component respectively. Here, the field size would be reduced by a factor of 8,

For x and y coarsening the coarser output cells represent the average value of all the calculation cells they represent. For z, the extra cells are simply skipped as only one calculation slice is present in memory at a time.

If the coarsening ratio doesn't cleanly divide the number of cells, some output cells at the edges will have to represent less calculation cells than the middle ones, causing the output to be shifted to the left by about half an output cell.

  • Small enough (< few 100s of lines), otherwise it should probably be split into smaller PRs
  • Tested (describe the tests in the PR description)
  • Runs on GPU (basic: the code compiles and run well with the new module)
  • Contains an automated test (checksum and/or comparison with theory)
  • Documented: all elements (classes and their members, functions, namespaces, etc.) are documented
  • Constified (All that can be const is const)
  • Code is clean (no unwanted comments, )
  • Style and code conventions are respected at the bottom of https://github.com/Hi-PACE/hipace
  • Proper label and GitHub project, if applicable

@AlexanderSinn AlexanderSinn added component: diagnostics About any types of diagnostics component: fields About 3D fields and slices, field solvers etc. labels Jul 30, 2021
Copy link
Member

@MaxThevenet MaxThevenet left a comment

Choose a reason for hiding this comment

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

Looks great, thanks for this PR!! See one inline comment below, on the interpolation method.

On another point: once we agree on the interpolation method, do you think you could add a test for these new diagnostics? For instance, for the blowout example, you could run:

  • 3D diagnostics with 2 4 8 coarsening
  • 3D diagnostics without coarsening, but do the coarsening in Python in post-processing.

and then check that the arrays are exactly the same. Then you could do the same with xz and yz slices.

Finally, do you think you could AMREX_ALWAYS_ASSERT_WITH_MESSAGE that the coarsening ratio divides the number of grid points in each direction? Better to restrict what the user does than to have unexpected behavior.

Comment on lines 280 to 291
amrex::Real field_value = 0;
int n_values = 0;
for (int j_c = 0; j_c < coarse_y && (j*coarse_y+j_c) < ncells_y; ++j_c) {
for (int i_c = 0; i_c < coarse_x && (i*coarse_x+i_c) < ncells_x; ++i_c) {
field_value += slice_array( i*coarse_x+i_c,
j*coarse_y+j_c,
k*coarse_z,
n+slice_comp);
++n_values;
}
}
full_array(i,j,k,n+full_comp) = field_value / std::max(n_values,1);
Copy link
Member

Choose a reason for hiding this comment

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

If I understand correctly, this is large-order smoothing interpolation (where the stencil length is coarse_* for direction *), where the same coefficient is applies to all point of the stencil. This causes significant smoothing of the data, which could be unwanted. Besides, we probably don't have to re-write interpolation functions, we already have some (originally written for particles, but also used for field interpolation in #572).

In this case, I would recommend to do sub-sampling rather than large-stencil interpolation, where we would ideally just take 1 point in each direction. Of course, this is impractical because we want to conserve the centering of data. For this I would recommend to do as shown on the image below: interpolate between the two nearest grid points (in each direction), which is the minimum to fix the centering (this is what we do for slice I/O). Do you think you could implement this instead? Of course this is valid for x and y, and for 3D slice I/O.

IMG_20210809_104848

@AlexanderSinn
Copy link
Member Author

I have now added subsampling for x, y and z as well as a test for coarse 3D IO similar to the existing test for slice IO. The subsampling can handle boundary values even if the coarsening ratio doesn’t cleanly divide the number of grid points. In that case the last cell(s) still in the box are taken. This enables an effective xy output by specifying xyz with a z-coarsening-ratio of >= to the number of grid points in z. The output slice would be z-coarsening-ratio/2.

Copy link
Member

@MaxThevenet MaxThevenet left a comment

Choose a reason for hiding this comment

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

That looks great to me, thanks a lot for this PR! Just one point: did you try to combine with with mesh refinement? Could you quickly test it? If it doesn't work as expected, that's not a problem, but we could assert somewhere that I/O coarsening and MR aren't used together, that would be a nice safeguard.

@MaxThevenet
Copy link
Member

In the meantime, I'll merge this PR, adding an assert for MR can be done in a subsequent PR. Thanks again!

@MaxThevenet MaxThevenet merged commit d435e90 into Hi-PACE:development Sep 7, 2021
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
component: diagnostics About any types of diagnostics component: fields About 3D fields and slices, field solvers etc.
Projects
None yet
Development

Successfully merging this pull request may close these issues.

2 participants