-
Notifications
You must be signed in to change notification settings - Fork 376
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
ReduceToPlaneMF #4384
Draft
WeiqunZhang
wants to merge
2
commits into
AMReX-Codes:development
Choose a base branch
from
WeiqunZhang:reduceplane
base: development
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.
Draft
ReduceToPlaneMF #4384
Conversation
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
Add a new Reduce to plane function that returns the results in a pair of MultiFabs. The first MultiFab is distributed the same as the original data MultiFab with each Box flatten to 1 cell. The Fabs of the first MultiFab contain box-local reduction results. The second MultiFab has only a single cell in the reduction direction, and it contains global reduction results distributed in the other two directions.
Test code #include <AMReX.H>
#include <AMReX_ParmParse.H>
#include <AMReX_Print.H>
#include <AMReX_MultiFabUtil.H>
using namespace amrex;
int main(int argc, char* argv[])
{
amrex::Initialize(argc, argv);
{
int ncomp = 1;
int n_cell = 64;
int max_grid_size = 32;
Box domain(IntVect(0), IntVect(n_cell-1));
BoxArray ba(domain);
ba.maxSize(max_grid_size);
MultiFab mf(ba, DistributionMapping{ba}, ncomp, 0);
mf.setVal(1);
auto const& ma = mf.const_arrays();
auto [mf2, mf2_unique] = ReduceToPlaneMF<MultiFab,ReduceOpSum,Real>
(2, domain, mf, [=] AMREX_GPU_DEVICE (int b, int i, int j, int k)
{
return ma[b](i,j,k);
});
printCell(mf2, IntVect(3,4,0));
printCell(mf2_unique, IntVect(3,4,0));
// auto phi = poisson_solver(mf2_unique);
// mf2.ParallelCopy(phi);
}
amrex::Finalize();
} |
6 tasks
WeiqunZhang
commented
Mar 20, 2025
ba2d.maxSize(ba[0].length()); | ||
DistributionMapping dm2d(ba2d); | ||
FA mf2d(ba2d, dm2d, 1, 0); | ||
mf2d.ParallelAdd(tmpfa); |
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.
This is only correct if it's reduce sum. We don't have a way for min and max yet. So maybe we can add a static assert there.
Sign up for free
to join this conversation on GitHub.
Already have an account?
Sign in to comment
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.
Add a new Reduce to plane function that returns the results in a pair of MultiFabs. The first MultiFab is distributed the same as the original data MultiFab with each Box flatten to 1 cell. The Fabs of the first MultiFab contain box-local reduction results. The second MultiFab has only a single cell in the reduction direction, and it contains global reduction results distributed in the other two directions.