Skip to content

Conversation

@WeiqunZhang
Copy link
Member

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.

@WeiqunZhang
Copy link
Member Author

WeiqunZhang commented Mar 20, 2025

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();
}

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.
@asalmgren asalmgren merged commit 5e3a7a3 into AMReX-Codes:development Apr 23, 2025
75 checks passed
@ax3l ax3l mentioned this pull request Jul 15, 2025
5 tasks
WeiqunZhang added a commit that referenced this pull request Jul 15, 2025
## Summary

Fix an assert in `reduce_to_plane`, making sure the same staggering is
used for the number-of-node-points assert.

## Additional background

Follow-up to #4384

## Checklist

The proposed changes:
- [x] fix a bug or incorrect behavior in AMReX
- [ ] add new capabilities to AMReX
- [ ] changes answers in the test suite to more than roundoff level
- [ ] are likely to significantly affect the results of downstream AMReX
users
- [ ] include documentation in the code and/or rst files, if appropriate

Co-authored-by: Weiqun Zhang <[email protected]>
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Projects

None yet

Development

Successfully merging this pull request may close these issues.

3 participants