![]() |
Block-Structured AMR Software Framework
|
MPI data movement operations for amrex::FabArray based containers. More...
Functions | |
| template<typename BUF = value_type> | |
| void | amrex::FabArray< FAB >::FillBoundary (bool cross=false) |
| Copy on intersection within a FabArray. | |
| template<typename BUF = value_type> | |
| void | amrex::FabArray< FAB >::FillBoundary (const Periodicity &period, bool cross=false) |
| template<typename BUF = value_type> | |
| void | amrex::FabArray< FAB >::FillBoundary (int scomp, int ncomp, bool cross=false) |
| template<typename BUF = value_type> | |
| void | amrex::FabArray< FAB >::FillBoundary (int scomp, int ncomp, const Periodicity &period, bool cross=false) |
| template<typename BUF = value_type> | |
| void | amrex::FabArray< FAB >::FillBoundary (int scomp, int ncomp, const IntVect &nghost, const Periodicity &period, bool cross=false) |
| void | amrex::FabArray< FAB >::FillBoundaryAndSync (const Periodicity &period=Periodicity::NonPeriodic()) |
| Fill ghost cells and synchronize nodal data. | |
| void | amrex::FabArray< FAB >::FillBoundaryAndSync (int scomp, int ncomp, const IntVect &nghost, const Periodicity &period) |
| Fill ghost cells and synchronize nodal data. | |
| void | amrex::FabArray< FAB >::OverrideSync (const Periodicity &period=Periodicity::NonPeriodic()) |
| Synchronize nodal data. | |
| void | amrex::FabArray< FAB >::OverrideSync (int scomp, int ncomp, const Periodicity &period) |
| Synchronize nodal data. | |
| void | amrex::FabArray< FAB >::SumBoundary (const Periodicity &period=Periodicity::NonPeriodic(), bool deterministic=false) |
| Sum values in overlapped cells. | |
| void | amrex::FabArray< FAB >::SumBoundary (int scomp, int ncomp, const Periodicity &period=Periodicity::NonPeriodic(), bool deterministic=false) |
| void | amrex::FabArray< FAB >::SumBoundary (int scomp, int ncomp, IntVect const &nghost, const Periodicity &period=Periodicity::NonPeriodic(), bool deterministic=false) |
| Sum values in overlapped cells. | |
| void | amrex::FabArray< FAB >::SumBoundary (int scomp, int ncomp, IntVect const &src_nghost, IntVect const &dst_nghost, const Periodicity &period=Periodicity::NonPeriodic(), bool deterministic=false) |
| Sum values in overlapped cells. | |
| void | amrex::FabArray< FAB >::EnforcePeriodicity (const Periodicity &period) |
| Fill ghost cells with values from their corresponding cells across periodic boundaries, regardless of whether the corresponding cells are valid. | |
| void | amrex::FabArray< FAB >::EnforcePeriodicity (int scomp, int ncomp, const Periodicity &period) |
| void | amrex::FabArray< FAB >::EnforcePeriodicity (int scomp, int ncomp, const IntVect &nghost, const Periodicity &period) |
| template<typename F = FAB, std::enable_if_t< IsBaseFab< F >::value, int > = 0> | |
| F::value_type | amrex::FabArray< FAB >::norminf (int comp, int ncomp, IntVect const &nghost, bool local=false, bool ignore_covered=false) const |
| Return infinity norm. | |
| template<typename IFAB , typename F = FAB, std::enable_if_t< IsBaseFab< F >::value, int > = 0> | |
| F::value_type | amrex::FabArray< FAB >::norminf (FabArray< IFAB > const &mask, int comp, int ncomp, IntVect const &nghost, bool local=false) const |
| Return infinity norm in masked region. | |
| Real | amrex::MultiFab::norminf (int comp=0, int nghost=0, bool local=false, bool ignore_covered=false) const |
| Real | amrex::MultiFab::norminf (const iMultiFab &mask, int comp=0, int nghost=0, bool local=false) const |
| Vector< Real > | amrex::MultiFab::norminf (const Vector< int > &comps, int nghost=0, bool local=false, bool ignore_covered=false) const |
| Real | amrex::MultiFab::norm1 (int comp=0, int ngrow=0, bool local=false) const |
Returns the L1 norm of component comp over the MultiFab. ngrow ghost cells are used. | |
| Vector< Real > | amrex::MultiFab::norm1 (const Vector< int > &comps, int ngrow=0, bool local=false) const |
| Returns the L1 norm of each component of "comps" over the MultiFab. ngrow ghost cells are used. | |
| Real | amrex::MultiFab::norm2 (int comp=0) const |
Returns the L2 norm of component comp over the MultiFab. No ghost cells are used. | |
| Real | amrex::MultiFab::norm2 (int comp, int numcomp) const |
Returns the L2 norm of numcomp components starting from comp component over the MultiFab. No ghost cells are used. | |
| Real | amrex::MultiFab::norm2 (int comp, const Periodicity &period) const |
Returns the L2 norm of component comp over the MultiFab. No ghost cells are used. This version has no double counting for nodal data. | |
| Vector< Real > | amrex::MultiFab::norm2 (const Vector< int > &comps) const |
| Returns the L2 norm of each component of "comps" over the MultiFab. No ghost cells are used. | |
| Real | amrex::MultiFab::sum (int comp=0, bool local=false) const |
| Returns the sum of component "comp" over the MultiFab – no ghost cells are included. | |
| Real | amrex::MultiFab::sum (Box const ®ion, int comp=0, bool local=false) const |
| Returns the sum of component "comp" in the given "region". – no ghost cells are included. | |
| Real | amrex::MultiFab::sum_unique (int comp=0, bool local=false, const Periodicity &period=Periodicity::NonPeriodic()) const |
Same as sum with local =false, but for non-cell-centered data, this only adds non-unique points that are owned by multiple boxes once. | |
| IntVect | amrex::MultiFab::minIndex (int comp, int nghost=0) const |
| IntVect | amrex::MultiFab::maxIndex (int comp, int nghost=0) const |
| static Real | amrex::MultiFab::Dot (const MultiFab &x, int xcomp, const MultiFab &y, int ycomp, int numcomp, int nghost, bool local=false) |
| Returns the dot product of two MultiFabs. | |
| static Real | amrex::MultiFab::Dot (const MultiFab &x, int xcomp, int numcomp, int nghost, bool local=false) |
| Returns the dot product of a MultiFab with itself. | |
| static Real | amrex::MultiFab::Dot (const iMultiFab &mask, const MultiFab &x, int xcomp, const MultiFab &y, int ycomp, int numcomp, int nghost, bool local=false) |
MPI data movement operations for amrex::FabArray based containers.
These operations may use MPI when AMReX is built with MPI support. Central interfaces include:
|
static |
|
static |
Returns the dot product of two MultiFabs.
|
static |
Returns the dot product of a MultiFab with itself.
| void amrex::FabArray< FAB >::EnforcePeriodicity | ( | const Periodicity & | period | ) |
Fill ghost cells with values from their corresponding cells across periodic boundaries, regardless of whether the corresponding cells are valid.
This differs from FillBoundary, which only fills from valid cells, and does not fill from ghost cells. The BoxArray is allowed to be overlapping.
| void amrex::FabArray< FAB >::EnforcePeriodicity | ( | int | scomp, |
| int | ncomp, | ||
| const IntVect & | nghost, | ||
| const Periodicity & | period | ||
| ) |
| void amrex::FabArray< FAB >::EnforcePeriodicity | ( | int | scomp, |
| int | ncomp, | ||
| const Periodicity & | period | ||
| ) |
| void amrex::FabArray< FAB >::FillBoundary | ( | bool | cross = false | ) |
Copy on intersection within a FabArray.
Data is copied from valid regions to intersecting regions of definition. The purpose is to fill in the boundary regions of each FAB in the FabArray. If cross=true, corner cells are not filled. If the length of periodic is provided, periodic boundaries are also filled. Note that FabArray itself does not contains any periodicity information. FillBoundary expects that its cell-centered version of its BoxArray is non-overlapping.
| void amrex::FabArray< FAB >::FillBoundary | ( | const Periodicity & | period, |
| bool | cross = false |
||
| ) |
| void amrex::FabArray< FAB >::FillBoundary | ( | int | scomp, |
| int | ncomp, | ||
| bool | cross = false |
||
| ) |
Same as FillBoundary(), but only copies ncomp components starting at scomp.
| void amrex::FabArray< FAB >::FillBoundary | ( | int | scomp, |
| int | ncomp, | ||
| const IntVect & | nghost, | ||
| const Periodicity & | period, | ||
| bool | cross = false |
||
| ) |
| void amrex::FabArray< FAB >::FillBoundary | ( | int | scomp, |
| int | ncomp, | ||
| const Periodicity & | period, | ||
| bool | cross = false |
||
| ) |
| void amrex::FabArray< FAB >::FillBoundaryAndSync | ( | const Periodicity & | period = Periodicity::NonPeriodic() | ) |
Fill ghost cells and synchronize nodal data.
Ghost regions are filled with data from the intersecting valid regions. The synchronization will override valid regions by the intersecting valid regions with a higher precedence. The smaller the global box index is, the higher precedence the box has. With periodic boundaries, for cells in the same box, those near the lower corner have higher precedence than those near the upper corner.
| period | periodic length if it's non-zero |
| void amrex::FabArray< FAB >::FillBoundaryAndSync | ( | int | scomp, |
| int | ncomp, | ||
| const IntVect & | nghost, | ||
| const Periodicity & | period | ||
| ) |
Fill ghost cells and synchronize nodal data.
Ghost regions are filled with data from the intersecting valid regions. The synchronization will override valid regions by the intersecting valid regions with a higher precedence. The smaller the global box index is, the higher precedence the box has. With periodic boundaries, for cells in the same box, those near the lower corner have higher precedence than those near the upper corner.
| scomp | starting component |
| ncomp | number of components |
| nghost | number of ghost cells to fill |
| period | periodic length if it's non-zero |
| Vector< Real > amrex::MultiFab::norm1 | ( | const Vector< int > & | comps, |
| int | ngrow = 0, |
||
| bool | local = false |
||
| ) | const |
Returns the L1 norm of each component of "comps" over the MultiFab. ngrow ghost cells are used.
Returns the L1 norm of component comp over the MultiFab. ngrow ghost cells are used.
Returns the L2 norm of each component of "comps" over the MultiFab. No ghost cells are used.
| Real amrex::MultiFab::norm2 | ( | int | comp, |
| const Periodicity & | period | ||
| ) | const |
Returns the L2 norm of component comp over the MultiFab. No ghost cells are used. This version has no double counting for nodal data.
Returns the L2 norm of numcomp components starting from comp component over the MultiFab. No ghost cells are used.
Returns the L2 norm of component comp over the MultiFab. No ghost cells are used.
|
inline |
|
inline |
| F::value_type amrex::FabArray< FAB >::norminf | ( | FabArray< IFAB > const & | mask, |
| int | comp, | ||
| int | ncomp, | ||
| IntVect const & | nghost, | ||
| bool | local = false |
||
| ) | const |
Return infinity norm in masked region.
| mask | only mask=true region is included |
| comp | starting component |
| ncomp | number of components |
| nghost | number of ghost cells |
| local | If true, MPI communication is skipped. |
| F::value_type amrex::FabArray< FAB >::norminf | ( | int | comp, |
| int | ncomp, | ||
| IntVect const & | nghost, | ||
| bool | local = false, |
||
| bool | ignore_covered = false |
||
| ) | const |
Return infinity norm.
| comp | starting component |
| ncomp | number of components |
| nghost | number of ghost cells |
| local | If true, MPI communication is skipped. |
| ignore_covered | ignore covered cells. Only relevant for cell-centered EB data. |
|
inline |
| void amrex::FabArray< FAB >::OverrideSync | ( | const Periodicity & | period = Periodicity::NonPeriodic() | ) |
Synchronize nodal data.
The synchronization will override valid regions by the intersecting valid regions with a higher precedence. The smaller the global box index is, the higher precedence the box has. With periodic boundaries, for cells in the same box, those near the lower corner have higher precedence than those near the upper corner.
Example: Suppose there is a valid point that is shared by 4 boxes. Some operation assigns values a, b, c, and d to the 4 points in the 4 boxes. We then call SumBoundary to add the 4 values and store the result in the 4 boxes. The value in box 0 might be the result of a+b+c+d, the value in box 1 might be the result of b+c+d+a, etc. The resulting data could be out of sync across boxes due to roundoff errors. OverrideSync synchronizes the data by overriding all values with one value (e.g., the value from box 0).
| period | periodic length if it's non-zero |
| void amrex::FabArray< FAB >::OverrideSync | ( | int | scomp, |
| int | ncomp, | ||
| const Periodicity & | period | ||
| ) |
Synchronize nodal data.
The synchronization will override valid regions by the intersecting valid regions with a higher precedence. The smaller the global box index is, the higher precedence the box has. With periodic boundaries, for cells in the same box, those near the lower corner have higher precedence than those near the upper corner.
| scomp | starting component |
| ncomp | number of components |
| period | periodic length if it's non-zero |
Returns the sum of component "comp" in the given "region". – no ghost cells are included.
Returns the sum of component "comp" over the MultiFab – no ghost cells are included.
| Real amrex::MultiFab::sum_unique | ( | int | comp = 0, |
| bool | local = false, |
||
| const Periodicity & | period = Periodicity::NonPeriodic() |
||
| ) | const |
Same as sum with local =false, but for non-cell-centered data, this only adds non-unique points that are owned by multiple boxes once.
| void amrex::FabArray< FAB >::SumBoundary | ( | const Periodicity & | period = Periodicity::NonPeriodic(), |
| bool | deterministic = false |
||
| ) |
Sum values in overlapped cells.
The destination is limited to valid cells. Note that being deterministic is not the same as data consistency on shared nodes. If you want data consistency on shared nodes, you can call OverrideSync after this. Sometimes, you can also wait till you need to do ghost cell exchanges, and you can then call FillBoundaryAndSync to do both.
| void amrex::FabArray< FAB >::SumBoundary | ( | int | scomp, |
| int | ncomp, | ||
| const Periodicity & | period = Periodicity::NonPeriodic(), |
||
| bool | deterministic = false |
||
| ) |
| void amrex::FabArray< FAB >::SumBoundary | ( | int | scomp, |
| int | ncomp, | ||
| IntVect const & | nghost, | ||
| const Periodicity & | period = Periodicity::NonPeriodic(), |
||
| bool | deterministic = false |
||
| ) |
Sum values in overlapped cells.
The destination is limited to valid + ngrow cells. Note that being deterministic is not the same as data consistency on shared nodes. If you want data consistency on shared nodes, you can call OverrideSync after this. Sometimes, you can also wait till you need to do ghost cell exchanges, and you can then call FillBoundaryAndSync to do both.
| void amrex::FabArray< FAB >::SumBoundary | ( | int | scomp, |
| int | ncomp, | ||
| IntVect const & | src_nghost, | ||
| IntVect const & | dst_nghost, | ||
| const Periodicity & | period = Periodicity::NonPeriodic(), |
||
| bool | deterministic = false |
||
| ) |
Sum values in overlapped cells.
For computing the overlap, the dst is grown by dst_ngrow, while the src uses src_ngrow. Note that being deterministic is not the same as data consistency on shared nodes. If you want data consistency on shared nodes, you can call OverrideSync after this. Sometimes, you can also wait till you need to do ghost cell exchanges, and you can then call FillBoundaryAndSync to do both.