Python API

Imports

pyAMReX provides the public imports amrex.space1d, amrex.space2d and amrex.space3d, mirroring the compile-time option AMReX_SPACEDIM.

Due to limitations in AMReX, currently, only one of the imports can be used at a time in the same Python process. For example:

import amrex.space3d as amr

A 1D or 2D AMReX run needs its own Python process. Another dimensionality cannot be imported into the same Python process after choosing a specific dimensionality for import.

For brevity, below the 3D APIs are shown. pyAMReX classes and functions follow the same structure as the C++ AMReX APIs.

Base

class amrex.space3d.AMReX[source]
static empty() bool[source]
static erase(arg0: AMReX) None[source]
static size() int[source]
static top() AMReX[source]
class amrex.space3d.Config[source]
amrex_version: ClassVar[str] = '26.06'
gpu_backend = None
have_eb: ClassVar[bool] = True
have_gpu: ClassVar[bool] = False
have_mpi: ClassVar[bool] = True
have_omp: ClassVar[bool] = False
have_simd: ClassVar[bool] = False
precision: ClassVar[str] = 'DOUBLE'
precision_particles: ClassVar[str] = 'DOUBLE'
simd_size: ClassVar[int] = 1
spacedim: ClassVar[int] = 3
verbose: ClassVar[int] = 1
amrex.space3d.initialize(arg0: list) AMReX[source]

Initialize AMReX library

amrex.space3d.finalize(*args, **kwds)

Helper for @overload to raise when called.

amrex.space3d.initialized() bool[source]

Returns true if there are any currently-active and initialized AMReX instances (i.e. one for which amrex::Initialize has been called, and amrex::Finalize has not). Otherwise false.

amrex.space3d.size() int[source]

The amr stack size, the number of amr instances pushed.

class amrex.space3d.Arena[source]
static finalize() None[source]
has_free_device_memory(sz: SupportsInt | SupportsIndex) bool[source]

Does the device have enough free memory for allocating this much memory? For CPU builds, this always return true.

static initialize(arg0: bool) None[source]
property is_device: bool
property is_device_accessible: bool
property is_host_accessible: bool
property is_managed: bool
property is_pinned: bool
static print_usage(arg0: bool) None[source]
static print_usage_to_files(filename: str, message: str) None[source]
class amrex.space3d.Direction[source]
class amrex.space3d.CoordSys[source]
class amrex.space3d.CoordSys(arg0: CoordSys)
Coord() CoordType[source]
CoordInt() int[source]
class CoordType(value)[source]
RZ: ClassVar[CoordType]
SPHERICAL: ClassVar[CoordType]
cartesian: ClassVar[CoordType]
undef: ClassVar[CoordType]
IsCartesian() bool[source]
IsRZ() bool[source]
IsSPHERICAL() bool[source]
RZ: ClassVar[CoordType]
SPHERICAL: ClassVar[CoordType]
SetCoord(arg0: CoordType) None[source]
cartesian: ClassVar[CoordType]
ok() bool[source]
undef: ClassVar[CoordType]
class amrex.space3d.DistributionMapping[source]
class amrex.space3d.DistributionMapping(arg0: DistributionMapping)
class amrex.space3d.DistributionMapping(arg0: Vector_int)
class amrex.space3d.DistributionMapping(boxes: BoxArray)
class amrex.space3d.DistributionMapping(boxes: BoxArray, nprocs: SupportsInt | SupportsIndex)
ProcessorMap() Vector_int[source]
property capacity: int
define(**kwds)

Helper for @overload to raise when called.

property empty: bool
property size: int
class amrex.space3d.GeometryData[source]
CellSize(**kwds)

Helper for @overload to raise when called.

Coord() int[source]

return integer coordinate type

Domain() Box[source]

Returns our rectangular domain

ProbHi(**kwds)

Helper for @overload to raise when called.

ProbLo(**kwds)

Helper for @overload to raise when called.

property coord: int

The Coordinates type.

property domain: Box

The index domain.

property dx: Annotated[list[float], 'FixedSize(3)']

The cellsize for each coordinate direction.

isPeriodic(**kwds)

Helper for @overload to raise when called.

property is_periodic: Annotated[list[int], 'FixedSize(3)']

Returns whether the domain is periodic in each coordinate direction.

property prob_domain: RealBox

The problem domain (real).

class amrex.space3d.Geometry[source]
class amrex.space3d.Geometry(dom: Box, rb: RealBox, coord: SupportsInt | SupportsIndex, is_per: Annotated[collections.abc.Sequence[SupportsInt | SupportsIndex], 'FixedSize(3)'])
ProbHi(**kwds)

Helper for @overload to raise when called.

ProbLength(arg0: SupportsInt | SupportsIndex) float[source]

length of problem domain in specified dimension

ProbLo(**kwds)

Helper for @overload to raise when called.

ProbSize() float[source]

the overall size of the domain

ResetDefaultCoord() None[source]

Reset default coord of Geometry class with an Array of int

ResetDefaultPeriodicity() None[source]

Reset default periodicity of Geometry class with an Array of int

ResetDefaultProbDomain() None[source]

Reset default problem domain of Geometry class with a RealBox

coarsen(rr: IntVect3D) None[source]
data() GeometryData[source]

Returns non-static copy of geometry’s stored data

define(dom: Box, rb: RealBox, coord: SupportsInt | SupportsIndex, is_per: Annotated[Sequence[SupportsInt | SupportsIndex], 'FixedSize(3)']) None[source]

Set geometry

property domain: Box

The rectangular domain (index space).

growNonPeriodicDomain(**kwds)

Helper for @overload to raise when called.

growPeriodicDomain(**kwds)

Helper for @overload to raise when called.

insideRoundOffDomain(x: SupportsFloat | SupportsIndex, y: SupportsFloat | SupportsIndex, z: SupportsFloat | SupportsIndex) bool[source]

Returns true if a point is inside the roundoff domain. All particles with positions inside the roundoff domain are sure to be mapped to cells inside the Domain() box. Note that the same need not be true for all points inside ProbDomain()

isAllPeriodic() bool[source]

Is domain periodic in all directions?

isAnyPeriodic() bool[source]

Is domain periodic in any direction?

isPeriodic(**kwds)

Helper for @overload to raise when called.

outsideRoundOffDomain(x: SupportsFloat | SupportsIndex, y: SupportsFloat | SupportsIndex, z: SupportsFloat | SupportsIndex) bool[source]

Returns true if a point is outside the roundoff domain. All particles with positions inside the roundoff domain are sure to be mapped to cells inside the Domain() box. Note that the same need not be true for all points inside ProbDomain()

period(dir: SupportsInt | SupportsIndex) int[source]

Return the period in the specified direction

periodicity(**kwds)

Helper for @overload to raise when called.

property prob_domain: RealBox

The problem domain (real).

refine(rr: IntVect3D) None[source]
setPeriodicity(period: Annotated[Sequence[SupportsInt | SupportsIndex], 'FixedSize(3)']) Annotated[list[int], 'FixedSize(3)'][source]

Set periodicity flags and return the old flags. Note that, unlike Periodicity class, the flags are just boolean.

amrex.space3d.ParallelDescriptor.IOProcessor() bool[source]
amrex.space3d.ParallelDescriptor.IOProcessorNumber() int[source]
amrex.space3d.ParallelDescriptor.MyProc() int[source]
amrex.space3d.ParallelDescriptor.NProcs() int[source]
class amrex.space3d.Periodicity[source]
class amrex.space3d.Periodicity(arg0: IntVect3D)
property domain: Box

Cell-centered domain Box “infinitely” long in non-periodic directions.

property is_all_periodic: bool
property is_any_periodic: bool
is_periodic(dir: SupportsInt | SupportsIndex) bool[source]
static non_periodic() Periodicity[source]

Return the Periodicity object that is not periodic in any direction

property shift_IntVect: list[IntVect3D]
class amrex.space3d.RealBox[source]
class amrex.space3d.RealBox(x_lo: SupportsFloat | SupportsIndex, y_lo: SupportsFloat | SupportsIndex, z_lo: SupportsFloat | SupportsIndex, x_hi: SupportsFloat | SupportsIndex, y_hi: SupportsFloat | SupportsIndex, z_hi: SupportsFloat | SupportsIndex)
class amrex.space3d.RealBox(a_lo: Annotated[collections.abc.Sequence[SupportsFloat | SupportsIndex], 'FixedSize(3)'], a_hi: Annotated[collections.abc.Sequence[SupportsFloat | SupportsIndex], 'FixedSize(3)'])
class amrex.space3d.RealBox(bx: Box, dx: Annotated[collections.abc.Sequence[SupportsFloat | SupportsIndex], 'FixedSize(3)'], base: Annotated[collections.abc.Sequence[SupportsFloat | SupportsIndex], 'FixedSize(3)'])
contains(**kwds)

Helper for @overload to raise when called.

hi(**kwds)

Helper for @overload to raise when called.

intersects(arg0: RealBox) bool[source]

determine if box intersects with a box

length(arg0: SupportsInt | SupportsIndex) float[source]
lo(**kwds)

Helper for @overload to raise when called.

ok() bool[source]

Determine if RealBox satisfies xlo[i]<xhi[i] for i=0,1,...,AMREX_SPACEDIM.

setHi(**kwds)

Helper for @overload to raise when called.

setLo(**kwds)

Helper for @overload to raise when called.

volume() float[source]
property xhi: Annotated[list[float], 'FixedSize(3)']
property xlo: Annotated[list[float], 'FixedSize(3)']
amrex.space3d.AlmostEqual(rb1: RealBox, rb2: RealBox, eps: SupportsFloat | SupportsIndex = 0.0) bool[source]

Determine if two boxes are equal to within a tolerance

Indexing: Box, IntVect and IndexType

Corresponding AMReX documentation.

amrex.space3d.IntVect

alias of IntVect3D

class amrex.space3d.Box(small: IntVect3D, big: IntVect3D)[source]
class amrex.space3d.Box(small: IntVect3D, big: IntVect3D, typ: IntVect3D)
class amrex.space3d.Box(small: IntVect3D, big: IntVect3D, t: IndexType)
class amrex.space3d.Box(small: Annotated[collections.abc.Sequence[SupportsInt | SupportsIndex], 'FixedSize(3)'], big: Annotated[collections.abc.Sequence[SupportsInt | SupportsIndex], 'FixedSize(3)'])
class amrex.space3d.Box(small: Annotated[collections.abc.Sequence[SupportsInt | SupportsIndex], 'FixedSize(3)'], big: Annotated[collections.abc.Sequence[SupportsInt | SupportsIndex], 'FixedSize(3)'], t: IndexType)
begin(box: Box) Dim3[source]
big_end: IntVect3D
property cell_centered: bool

Returns true if Box is cell-centered in all indexing directions.

contains(p: IntVect3D) bool[source]

Returns true if argument is contained within Box.

convert(**kwds)

Helper for @overload to raise when called.

property d_num_pts: float
enclosed_cells(**kwds)

Helper for @overload to raise when called.

end(box: Box) Dim3[source]
grow(**kwds)

Helper for @overload to raise when called.

grow_high(**kwds)

Helper for @overload to raise when called.

grow_low(**kwds)

Helper for @overload to raise when called.

hi_vect: IntVect3D
intersects(b: Box) bool[source]

Returns true if Boxes have non-null intersections. It is an error if the Boxes have different types.

property is_empty: bool
property is_square: bool
property ix_type: IndexType
lbound(arg0: Box) Dim3[source]
length(**kwds)

Helper for @overload to raise when called.

lo_vect: IntVect3D
make_slab(direction: SupportsInt | SupportsIndex, slab_index: SupportsInt | SupportsIndex) Box[source]

Flatten the box in one direction.

normalize() None[source]
numPts() int[source]

Return the number of points in the Box.

property num_pts: int
property ok: bool
same_size(b: Box) bool[source]

Returns true is Boxes same size, ie translates of each other,. It is an error if they have different types.

same_type(b: Box) bool[source]

Returns true if Boxes have same type.

shift(**kwds)

Helper for @overload to raise when called.

property size: IntVect3D
small_end: IntVect3D
strictly_contains(p: IntVect3D) bool[source]

Returns true if argument is strictly contained within Box.

surrounding_nodes(**kwds)

Helper for @overload to raise when called.

property the_unit_box: Box
property type: IntVect3D
ubound(arg0: Box) Dim3[source]
property volume: int
class amrex.space3d.BoxArray[source]
class amrex.space3d.BoxArray(arg0: BoxArray)
class amrex.space3d.BoxArray(arg0: Box)
class amrex.space3d.BoxArray(arg0: Vector_Box)
property capacity: int
cell_equal(arg0: BoxArray) bool[source]
clear() None[source]
coarsen(**kwds)

Helper for @overload to raise when called.

coarsenable(**kwds)

Helper for @overload to raise when called.

convert(**kwds)

Helper for @overload to raise when called.

property d_numPts: float
define(arg0: Box) None[source]
property empty: bool
enclosed_cells(**kwds)

Helper for @overload to raise when called.

get(arg0: SupportsInt | SupportsIndex) Box[source]
ix_type() IndexType[source]
max_size(**kwds)

Helper for @overload to raise when called.

minimal_box() Box[source]
property numPts: int
refine(**kwds)

Helper for @overload to raise when called.

resize(arg0: SupportsInt | SupportsIndex) None[source]
property size: int
surroundingNodes(**kwds)

Helper for @overload to raise when called.

class amrex.space3d.Dim3(arg0: SupportsInt | SupportsIndex, arg1: SupportsInt | SupportsIndex, arg2: SupportsInt | SupportsIndex)[source]
property x: int
property y: int
property z: int
class amrex.space3d.XDim3(arg0: SupportsFloat | SupportsIndex, arg1: SupportsFloat | SupportsIndex, arg2: SupportsFloat | SupportsIndex)[source]
property x: float
property y: float
property z: float
class amrex.space3d.IndexType[source]
class amrex.space3d.IndexType(arg0: IndexType)
class amrex.space3d.IndexType(arg0: IndexType.CellIndex, arg1: IndexType.CellIndex, arg2: IndexType.CellIndex)
CELL: ClassVar[CellIndex]
class CellIndex(value)[source]
CELL: ClassVar[CellIndex]
NODE: ClassVar[CellIndex]
NODE: ClassVar[CellIndex]
any() bool[source]
cell_centered(**kwds)

Helper for @overload to raise when called.

static cell_type() IndexType[source]
clear() None[source]
flip(arg0: SupportsInt | SupportsIndex) None[source]
ix_type(**kwds)

Helper for @overload to raise when called.

node_centered(**kwds)

Helper for @overload to raise when called.

static node_type() IndexType[source]
ok() bool[source]
set(arg0: SupportsInt | SupportsIndex) None[source]
set_type(arg0: SupportsInt | SupportsIndex, arg1: CellIndex) None[source]
setall() None[source]
test(arg0: SupportsInt | SupportsIndex) bool[source]
to_IntVect() IntVect3D[source]
unset(arg0: SupportsInt | SupportsIndex) None[source]

Vectors

class amrex.space3d.RealVect[source]
class amrex.space3d.RealVect(arg0: SupportsFloat | SupportsIndex, arg1: SupportsFloat | SupportsIndex, arg2: SupportsFloat | SupportsIndex)
class amrex.space3d.RealVect(arg0: IntVect3D)
class amrex.space3d.RealVect(arg0: collections.abc.Sequence[SupportsFloat | SupportsIndex])
class amrex.space3d.RealVect(arg0: SupportsFloat | SupportsIndex)
BASISREALV() RealVect[source]

return basis vector in given coordinate direction

ceil() IntVect3D[source]

Return an IntVect whose components are the std::ceil of the vector components

crossProduct(arg0: RealVect) RealVect[source]

Return cross product of this vector with another

dotProduct(arg0: RealVect) float[source]

Return dot product of this vector with another

floor() IntVect3D[source]

Return an IntVect whose components are the std::floor of the vector components

max(arg0: RealVect) RealVect[source]

Replace vector with the component-wise maxima of this vector and another

maxDir(arg0: bool) int[source]

direction or index of maximum value of this vector

min(arg0: RealVect) RealVect[source]

Replace vector with the component-wise minima of this vector and another

minDir(arg0: bool) int[source]

direction or index of minimum value of this vector

property product: float

Product of entries of this vector

property radSquared: float

Length squared of this vector

round() IntVect3D[source]

Return an IntVect whose components are the std::round of the vector components

scale(arg0: SupportsFloat | SupportsIndex) RealVect[source]

Multiplify each component of this vector by a scalar

property sum: float

Sum of the components of this vector

static unit_vector() RealVect[source]
property vectorLength: float

Length or 2-Norm of this vector

static zero_vector() RealVect[source]
amrex.space3d.min(arg0: RealVect, arg1: RealVect) RealVect[source]
amrex.space3d.max(arg0: RealVect, arg1: RealVect) RealVect[source]

amrex::Vector<T> is implemented for many types, e.g.,

class amrex.space3d.Vector_Real[source]
class amrex.space3d.Vector_Real(arg0: Vector_Real)
class amrex.space3d.Vector_Real(arg0: Iterable)
class amrex.space3d.Vector_Real
class amrex.space3d.Vector_Real(arg0: Vector_Real)
append(x: SupportsFloat | SupportsIndex) None[source]

Add an item to the end of the list

clear() None[source]

Clear the contents

count(x: SupportsFloat | SupportsIndex) int[source]

Return the number of times x appears in the list

extend(**kwds)

Helper for @overload to raise when called.

insert(i: SupportsInt | SupportsIndex, x: SupportsFloat | SupportsIndex) None[source]

Insert an item at a given position.

pop(**kwds)

Helper for @overload to raise when called.

remove(x: SupportsFloat | SupportsIndex) None[source]

Remove the first item from the list whose value is x. It is an error if there is no such item.

size() int[source]
class amrex.space3d.Vector_int[source]
class amrex.space3d.Vector_int(arg0: Vector_int)
class amrex.space3d.Vector_int(arg0: Iterable)
class amrex.space3d.Vector_int
class amrex.space3d.Vector_int(arg0: Vector_int)
append(x: SupportsInt | SupportsIndex) None[source]

Add an item to the end of the list

clear() None[source]

Clear the contents

count(x: SupportsInt | SupportsIndex) int[source]

Return the number of times x appears in the list

extend(**kwds)

Helper for @overload to raise when called.

insert(i: SupportsInt | SupportsIndex, x: SupportsInt | SupportsIndex) None[source]

Insert an item at a given position.

pop(**kwds)

Helper for @overload to raise when called.

remove(x: SupportsInt | SupportsIndex) None[source]

Remove the first item from the list whose value is x. It is an error if there is no such item.

size() int[source]
class amrex.space3d.Vector_string[source]
class amrex.space3d.Vector_string(arg0: Vector_string)
class amrex.space3d.Vector_string(arg0: Iterable)
class amrex.space3d.Vector_string
class amrex.space3d.Vector_string(arg0: Vector_string)
append(x: str) None[source]

Add an item to the end of the list

clear() None[source]

Clear the contents

count(x: str) int[source]

Return the number of times x appears in the list

extend(**kwds)

Helper for @overload to raise when called.

insert(i: SupportsInt | SupportsIndex, x: str) None[source]

Insert an item at a given position.

pop(**kwds)

Helper for @overload to raise when called.

remove(x: str) None[source]

Remove the first item from the list whose value is x. It is an error if there is no such item.

size() int[source]

Data Containers

amrex::Array4<T> is implemented for many floating point and integer types, e.g.,

class amrex.space3d.Array4_double[source]
class amrex.space3d.Array4_double(arg0: Array4_double)
class amrex.space3d.Array4_double(arg0: Array4_double, arg1: SupportsInt | SupportsIndex)
class amrex.space3d.Array4_double(arg0: Array4_double, arg1: SupportsInt | SupportsIndex, arg2: SupportsInt | SupportsIndex)
class amrex.space3d.Array4_double(arg0: Annotated[numpy.typing.ArrayLike, numpy.float64])
contains(**kwds)

Helper for @overload to raise when called.

property nComp: int
property num_comp: int
property size: int
to_cupy(copy=False, order='F')[source]

Provide a CuPy view into an Array4.

This includes ngrow guard cells of the box.

Note on the order of indices: By default, this is as in AMReX in Fortran contiguous order, indexing as x,y,z. This has performance implications for use in external libraries such as cupy. The order=”C” option will index as z,y,x and perform better with cupy. https://github.com/AMReX-Codes/pyamrex/issues/55#issuecomment-1579610074

Parameters

selfamrex.Array4_*

An Array4 class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

orderstring, optional

F order (default) or C. C is faster with external libraries.

Returns

cupy.array

A cupy n-dimensional array.

Raises

ImportError

Raises an exception if cupy is not installed

to_host() ndarray[tuple[Any, ...], dtype[float64]][source]
to_numpy(copy=False, order='F')[source]

Provide a NumPy view into an Array4.

This includes ngrow guard cells of the box.

Note on the order of indices: By default, this is as in AMReX in Fortran contiguous order, indexing as x,y,z. This has performance implications for use in external libraries such as cupy. The order=”C” option will index as z,y,x and perform better with cupy. https://github.com/AMReX-Codes/pyamrex/issues/55#issuecomment-1579610074

Parameters

selfamrex.Array4_*

An Array4 class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

orderstring, optional

F order (default) or C. C is faster with external libraries.

Returns

np.array

A NumPy n-dimensional array.

to_xp(copy=False, order='F')[source]

Provide a NumPy or CuPy view into an Array4, depending on amr.Config.have_gpu .

This function is similar to CuPy’s xp naming suggestion for CPU/GPU agnostic code: https://docs.cupy.dev/en/stable/user_guide/basic.html#how-to-write-cpu-gpu-agnostic-code

This includes ngrow guard cells of the box.

Note on the order of indices: By default, this is as in AMReX in Fortran contiguous order, indexing as x,y,z. This has performance implications for use in external libraries such as cupy. The order=”C” option will index as z,y,x and perform better with cupy. https://github.com/AMReX-Codes/pyamrex/issues/55#issuecomment-1579610074

Parameters

selfamrex.Array4_*

An Array4 class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

orderstring, optional

F order (default) or C. C is faster with external libraries.

Returns

xp.array

A NumPy or CuPy n-dimensional array.

class amrex.space3d.BaseFab_Real[source]
class amrex.space3d.BaseFab_Real(arg0: Arena)
class amrex.space3d.BaseFab_Real(arg0: Box, arg1: SupportsInt | SupportsIndex, arg2: Arena)
class amrex.space3d.BaseFab_Real(arg0: Box, arg1: SupportsInt | SupportsIndex, arg2: SupportsFloat | SupportsIndex)
class amrex.space3d.BaseFab_Real(arg0: Box, arg1: SupportsInt | SupportsIndex, arg2: SupportsFloat | SupportsIndex)
class amrex.space3d.BaseFab_Real(arg0: Array4_double)
class amrex.space3d.BaseFab_Real(arg0: Array4_double, arg1: IndexType)
class amrex.space3d.BaseFab_Real(arg0: Array4_double_const)
class amrex.space3d.BaseFab_Real(arg0: Array4_double_const, arg1: IndexType)
array() Array4_double[source]
big_end() IntVect3D[source]
box() Box[source]
clear() None[source]
const_array() Array4_double_const[source]
hi_vect() int[source]
is_allocated() bool[source]
length() IntVect3D[source]
lo_vect() int[source]
n_bytes(**kwds)

Helper for @overload to raise when called.

n_bytes_owned() int[source]
n_comp() int[source]
num_pts() int[source]
resize(arg0: Box, arg1: SupportsInt | SupportsIndex, arg2: Arena) None[source]
size() int[source]
small_end() IntVect3D[source]
to_host() BaseFab_Real[source]
class amrex.space3d.FArrayBox[source]
class amrex.space3d.FArrayBox(arg0: Arena)
class amrex.space3d.FArrayBox(arg0: Box, arg1: SupportsInt | SupportsIndex, arg2: Arena)
class amrex.space3d.FArrayBox(arg0: Box, arg1: SupportsInt | SupportsIndex, arg2: bool, arg3: bool, arg4: Arena)
class amrex.space3d.FArrayBox(arg0: Box, arg1: SupportsInt | SupportsIndex, arg2: SupportsFloat | SupportsIndex)
class amrex.space3d.FArrayBox(arg0: Box, arg1: SupportsInt | SupportsIndex, arg2: SupportsFloat | SupportsIndex)
class amrex.space3d.FArrayBox(arg0: Array4_double)
class amrex.space3d.FArrayBox(arg0: Array4_double, arg1: IndexType)
class amrex.space3d.FArrayBox(arg0: Array4_double_const)
class amrex.space3d.FArrayBox(arg0: Array4_double_const, arg1: IndexType)
class amrex.space3d.MultiFab[source]
class amrex.space3d.MultiFab(a: Arena)
class amrex.space3d.MultiFab(bxs: BoxArray, dm: DistributionMapping, ncomp: SupportsInt | SupportsIndex, ngrow: SupportsInt | SupportsIndex, info: MFInfo, factory: FabFactory_FArrayBox)
class amrex.space3d.MultiFab(bxs: BoxArray, dm: DistributionMapping, ncomp: SupportsInt | SupportsIndex, ngrow: SupportsInt | SupportsIndex, info: MFInfo)
class amrex.space3d.MultiFab(bxs: BoxArray, dm: DistributionMapping, ncomp: SupportsInt | SupportsIndex, ngrow: SupportsInt | SupportsIndex)
class amrex.space3d.MultiFab(bxs: BoxArray, dm: DistributionMapping, ncomp: SupportsInt | SupportsIndex, ngrow: IntVect3D, info: MFInfo)
class amrex.space3d.MultiFab(bxs: BoxArray, dm: DistributionMapping, ncomp: SupportsInt | SupportsIndex, ngrow: IntVect3D, info: MFInfo, factory: FabFactory_FArrayBox)
class amrex.space3d.MultiFab(bxs: BoxArray, dm: DistributionMapping, ncomp: SupportsInt | SupportsIndex, ngrow: IntVect3D)
add(**kwds)

Helper for @overload to raise when called.

add_product(**kwds)

Helper for @overload to raise when called.

average_sync(arg0: Periodicity) None[source]
box_array() BoxArray[source]
contains_inf(**kwds)

Helper for @overload to raise when called.

contains_nan(**kwds)

Helper for @overload to raise when called.

copy()[source]

Create a copy of this MultiFab, using the same Arena.

Parameters

selfamrex.MultiFab

A MultiFab class in pyAMReX

Returns

amrex.MultiFab

A copy of this MultiFab.

copymf(**kwds)

Helper for @overload to raise when called.

divi(mf: MultiFab, strt_comp: SupportsInt | SupportsIndex, num_comp: SupportsInt | SupportsIndex, nghost: SupportsInt | SupportsIndex = 0) None[source]

This function divides the values of the cells in mf from the corresponding cells of this MultiFab. mf is required to have the same BoxArray or “valid region” as this MultiFab. The division is done only to num_comp components, starting with component number strt_comp. The parameter nghost specifies the number of boundary cells that will be modified. If nghost == 0, only the valid region of each FArrayBox will be modified. Note, nothing is done to protect against divide by zero.

divide(**kwds)

Helper for @overload to raise when called.

dm() DistributionMapping[source]
dot(**kwds)

Helper for @overload to raise when called.

static finalize() None[source]
imesh(idir, include_ghosts=False)[source]
Returns the integer mesh along the specified direction with the appropriate centering.

This is the location of the data points in grid cell units.

selfamrex.MultiFab

A MultiFab class in pyAMReX

directioninteger

Zero based direction number. In a typical Cartesian case, 0 would be ‘x’ direction.

include_ghostsbool, default=False

Whether or not ghost cells are included in the mesh.

static initialize() None[source]
invert(**kwds)

Helper for @overload to raise when called.

lin_comb(a: SupportsFloat | SupportsIndex, x: MultiFab, x_comp: SupportsInt | SupportsIndex, b: SupportsFloat | SupportsIndex, y: MultiFab, y_comp: SupportsInt | SupportsIndex, comp: SupportsInt | SupportsIndex, numcomp: SupportsInt | SupportsIndex, nghost: SupportsInt | SupportsIndex) None[source]

self = a * x + b * y

max(**kwds)

Helper for @overload to raise when called.

maxIndex(arg0: SupportsInt | SupportsIndex, arg1: SupportsInt | SupportsIndex) IntVect3D[source]
min(**kwds)

Helper for @overload to raise when called.

minIndex(arg0: SupportsInt | SupportsIndex, arg1: SupportsInt | SupportsIndex) IntVect3D[source]
minus(mf: MultiFab, strt_comp: SupportsInt | SupportsIndex, num_comp: SupportsInt | SupportsIndex, nghost: SupportsInt | SupportsIndex = 0) None[source]

This function subtracts the values of the cells in mf from the corresponding cells of this MultiFab. mf is required to have the same BoxArray or “valid region” as this MultiFab. The subtraction is done only to num_comp components, starting with component number strt_comp. The parameter nghost specifies the number of boundary cells that will be modified. If nghost == 0, only the valid region of each FArrayBox will be modified.

mult(**kwds)

Helper for @overload to raise when called.

multiply(**kwds)

Helper for @overload to raise when called.

property n_comp: int
property n_grow_vect: IntVect3D

Return the grow factor (per direction) that defines the region of definition.

negate(**kwds)

Helper for @overload to raise when called.

norm0(**kwds)

Helper for @overload to raise when called.

norm1(**kwds)

Helper for @overload to raise when called.

norm2(**kwds)

Helper for @overload to raise when called.

norminf(arg0: SupportsInt | SupportsIndex, arg1: SupportsInt | SupportsIndex, arg2: bool, arg3: bool) float[source]
override_sync(arg0: iMultiFab, arg1: Periodicity) None[source]

Helper for @overload to raise when called.

plus(**kwds)

Helper for @overload to raise when called.

saxpy(a: SupportsFloat | SupportsIndex, src: MultiFab, srccomp: SupportsInt | SupportsIndex, comp: SupportsInt | SupportsIndex, numcomp: SupportsInt | SupportsIndex, nghost: SupportsInt | SupportsIndex) None[source]

self += a * src

property shape

Returns the shape of the global array

selfamrex.MultiFab

A MultiFab class in pyAMReX

include_ghostsbool, default=False

Whether or not ghost cells are included

property shape_with_ghosts

Returns the shape of the global array including ghost cells

selfamrex.MultiFab

A MultiFab class in pyAMReX

subtract(**kwds)

Helper for @overload to raise when called.

sum(**kwds)

Helper for @overload to raise when called.

sum_unique(**kwds)

Helper for @overload to raise when called.

swap(**kwds)

Helper for @overload to raise when called.

to_cupy(copy=False, order='F')[source]

Provide a CuPy view into a MultiFab.

This includes ngrow guard cells of each box.

Note on the order of indices: By default, this is as in AMReX in Fortran contiguous order, indexing as x,y,z. This has performance implications for use in external libraries such as cupy. The order=”C” option will index as z,y,x and perform better with cupy. https://github.com/AMReX-Codes/pyamrex/issues/55#issuecomment-1579610074

Parameters

selfamrex.MultiFab

A MultiFab class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

orderstring, optional

F order (default) or C. C is faster with external libraries.

Returns

list of cupy.array

A list of CuPy n-dimensional arrays, for each local block in the MultiFab.

Raises

ImportError

Raises an exception if cupy is not installed

to_numpy(copy=False, order='F')[source]

Provide a NumPy view into a MultiFab.

This includes ngrow guard cells of each box.

Note on the order of indices: By default, this is as in AMReX in Fortran contiguous order, indexing as x,y,z. This has performance implications for use in external libraries such as cupy. The order=”C” option will index as z,y,x and perform better with cupy. https://github.com/AMReX-Codes/pyamrex/issues/55#issuecomment-1579610074

Parameters

selfamrex.MultiFab

A MultiFab class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

orderstring, optional

F order (default) or C. C is faster with external libraries.

Returns

list of numpy.array

A list of NumPy n-dimensional arrays, for each local block in the MultiFab.

to_xp(copy=False, order='F')[source]

Provide a NumPy or CuPy view into a MultiFab, depending on amr.Config.have_gpu .

This function is similar to CuPy’s xp naming suggestion for CPU/GPU agnostic code: https://docs.cupy.dev/en/stable/user_guide/basic.html#how-to-write-cpu-gpu-agnostic-code

This includes ngrow guard cells of each box.

Note on the order of indices: By default, this is as in AMReX in Fortran contiguous order, indexing as x,y,z. This has performance implications for use in external libraries such as cupy. The order=”C” option will index as z,y,x and perform better with cupy. https://github.com/AMReX-Codes/pyamrex/issues/55#issuecomment-1579610074

Parameters

selfamrex.MultiFab

A MultiFab class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

orderstring, optional

F order (default) or C. C is faster with external libraries.

Returns

list of xp.array

A list of NumPy or CuPy n-dimensional arrays, for each local block in the MultiFab.

weighted_sync(arg0: MultiFab, arg1: Periodicity) None[source]
xpay(a: SupportsFloat | SupportsIndex, src: MultiFab, srccomp: SupportsInt | SupportsIndex, comp: SupportsInt | SupportsIndex, numcomp: SupportsInt | SupportsIndex, nghost: SupportsInt | SupportsIndex) None[source]

self = src + a * self

class amrex.space3d.MFInfo[source]
alloc: bool
arena: Arena
set_alloc(arg0: bool) MFInfo[source]
set_arena(arg0: Arena) MFInfo[source]
set_tag(arg0: str) None[source]
tags: Vector_string
class amrex.space3d.MFIter(arg0: FabArrayBase)[source]
class amrex.space3d.MFIter(arg0: FabArrayBase, arg1: MFItInfo)
class amrex.space3d.MFIter(arg0: MultiFab)
class amrex.space3d.MFIter(arg0: MultiFab, arg1: MFItInfo)
class amrex.space3d.MFIter(arg0: iMultiFab)
class amrex.space3d.MFIter(arg0: iMultiFab, arg1: MFItInfo)
fabbox() Box[source]
finalize() None[source]
grownnodaltilebox(**kwds)

Helper for @overload to raise when called.

growntilebox(ng: IntVect3D = -1000000) Box[source]
property index: int
property is_valid: bool
property length: int
nodaltilebox(dir: SupportsInt | SupportsIndex = -1) Box[source]
tilebox(**kwds)

Helper for @overload to raise when called.

validbox() Box[source]

amrex::PODVector<T, Allocator> is implemented for many allocators, e.g.,

class amrex.space3d.PODVector_real_arena[source]
class amrex.space3d.PODVector_real_arena(size: SupportsInt | SupportsIndex)
class amrex.space3d.PODVector_real_arena(other: PODVector_real_arena)

A plain-old-data (POD) vector of ‘real’ elements with ‘arena’ allocation.

assign(value: SupportsFloat | SupportsIndex) None[source]

assign the same value to every element

capacity() int[source]
clear() None[source]
empty() bool[source]
classmethod from_cupy(arr)[source]

Create a new PODVector from a CuPy array (or array-like).

Always copies the data into a newly allocated PODVector. Works for every allocator type: for host-only allocators the data is staged to the host through NumPy automatically.

Parameters

clstype

The PODVector type to construct.

arrarray_like

Input data, convertible to a CuPy array.

Returns

PODVector

A new PODVector with a copy of the data.

static from_numpy(arr: Annotated[_Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str], float64]) PODVector_real_arena[source]

Create a new PODVector from a NumPy array (or array-like).

Always copies the data into a newly allocated PODVector. The input is cast to the vector’s element type and made contiguous as needed. The copy into device-only memory uses an AMReX host-to-device copy and does not require CuPy.

Parameters

arrarray_like

Input data, convertible to a NumPy array.

Returns

PODVector

A new PODVector with a copy of the data.

classmethod from_xp(arr)[source]

Create a new PODVector from a NumPy or CuPy array, depending on amr.Config.have_gpu .

Always copies the data into a newly allocated PODVector. Unlike :meth:`to_xp`, a zero-copy view is not possible here because PODVector always owns its memory through its allocator.

This function is similar to CuPy’s xp naming suggestion for CPU/GPU agnostic code: https://docs.cupy.dev/en/stable/user_guide/basic.html#how-to-write-cpu-gpu-agnostic-code

Parameters

clstype

The PODVector type to construct.

arrarray_like

Input data (NumPy or CuPy array).

Returns

PODVector

A new PODVector with a copy of the data.

pop_back() None[source]
push_back(arg0: SupportsFloat | SupportsIndex) None[source]
reserve(capacity: SupportsInt | SupportsIndex, strategy: GrowthStrategy = 'GrowthStrategy.Poisson') None[source]
resize(**kwds)

Helper for @overload to raise when called.

shrink_to_fit() None[source]
size() int[source]
to_cupy(copy=False)[source]

Provide a CuPy view into a PODVector (e.g., RealVector, IntVector).

Parameters

selfamrex.PODVector_*

A PODVector class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

Returns

cupy.array

A 1D cupy array.

Raises

ImportError

Raises an exception if cupy is not installed

to_device() PODVector_real_std[source]

Copy this vector into a new amrex Gpu::DeviceVector (the arena allocator on GPU, std on CPU), transferring across memory spaces as needed. Mirrors to_host().

to_host() PODVector_real_pinned[source]

Copy this vector into a new pinned (host) PODVector. Mirrors to_device().

to_numpy(copy=False)[source]

Provide a NumPy view into a PODVector (e.g., RealVector, IntVector).

Parameters

selfamrex.PODVector_*

A PODVector class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

Returns

np.array

A 1D NumPy array.

to_xp(copy=False)[source]

Provide a NumPy or CuPy view into a PODVector (e.g., RealVector, IntVector), depending on amr.Config.have_gpu .

This function is similar to CuPy’s xp naming suggestion for CPU/GPU agnostic code: https://docs.cupy.dev/en/stable/user_guide/basic.html#how-to-write-cpu-gpu-agnostic-code

Parameters

selfamrex.PODVector_*

A PODVector class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

Returns

xp.array

A 1D NumPy or CuPy array.

class amrex.space3d.PODVector_int_pinned[source]
class amrex.space3d.PODVector_int_pinned(size: SupportsInt | SupportsIndex)
class amrex.space3d.PODVector_int_pinned(other: PODVector_int_pinned)

A plain-old-data (POD) vector of ‘int’ elements with ‘pinned’ allocation.

assign(value: SupportsInt | SupportsIndex) None[source]

assign the same value to every element

capacity() int[source]
clear() None[source]
empty() bool[source]
classmethod from_cupy(arr)[source]

Create a new PODVector from a CuPy array (or array-like).

Always copies the data into a newly allocated PODVector. Works for every allocator type: for host-only allocators the data is staged to the host through NumPy automatically.

Parameters

clstype

The PODVector type to construct.

arrarray_like

Input data, convertible to a CuPy array.

Returns

PODVector

A new PODVector with a copy of the data.

static from_numpy(arr: Annotated[_Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str], int32]) PODVector_int_pinned[source]

Create a new PODVector from a NumPy array (or array-like).

Always copies the data into a newly allocated PODVector. The input is cast to the vector’s element type and made contiguous as needed. The copy into device-only memory uses an AMReX host-to-device copy and does not require CuPy.

Parameters

arrarray_like

Input data, convertible to a NumPy array.

Returns

PODVector

A new PODVector with a copy of the data.

classmethod from_xp(arr)[source]

Create a new PODVector from a NumPy or CuPy array, depending on amr.Config.have_gpu .

Always copies the data into a newly allocated PODVector. Unlike :meth:`to_xp`, a zero-copy view is not possible here because PODVector always owns its memory through its allocator.

This function is similar to CuPy’s xp naming suggestion for CPU/GPU agnostic code: https://docs.cupy.dev/en/stable/user_guide/basic.html#how-to-write-cpu-gpu-agnostic-code

Parameters

clstype

The PODVector type to construct.

arrarray_like

Input data (NumPy or CuPy array).

Returns

PODVector

A new PODVector with a copy of the data.

pop_back() None[source]
push_back(arg0: SupportsInt | SupportsIndex) None[source]
reserve(capacity: SupportsInt | SupportsIndex, strategy: GrowthStrategy = 'GrowthStrategy.Poisson') None[source]
resize(**kwds)

Helper for @overload to raise when called.

shrink_to_fit() None[source]
size() int[source]
to_cupy(copy=False)[source]

Provide a CuPy view into a PODVector (e.g., RealVector, IntVector).

Parameters

selfamrex.PODVector_*

A PODVector class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

Returns

cupy.array

A 1D cupy array.

Raises

ImportError

Raises an exception if cupy is not installed

to_device() PODVector_int_std[source]

Copy this vector into a new amrex Gpu::DeviceVector (the arena allocator on GPU, std on CPU), transferring across memory spaces as needed. Mirrors to_host().

to_host() PODVector_int_pinned[source]

Copy this vector into a new pinned (host) PODVector. Mirrors to_device().

to_numpy(copy=False)[source]

Provide a NumPy view into a PODVector (e.g., RealVector, IntVector).

Parameters

selfamrex.PODVector_*

A PODVector class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

Returns

np.array

A 1D NumPy array.

to_xp(copy=False)[source]

Provide a NumPy or CuPy view into a PODVector (e.g., RealVector, IntVector), depending on amr.Config.have_gpu .

This function is similar to CuPy’s xp naming suggestion for CPU/GPU agnostic code: https://docs.cupy.dev/en/stable/user_guide/basic.html#how-to-write-cpu-gpu-agnostic-code

Parameters

selfamrex.PODVector_*

A PODVector class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

Returns

xp.array

A 1D NumPy or CuPy array.

Small Matrices and Vectors

class amrex.space3d.SmallMatrix_6x6_F_SI1_double[source]
class amrex.space3d.SmallMatrix_6x6_F_SI1_double(arg0: SmallMatrix_6x6_F_SI1_double)
class amrex.space3d.SmallMatrix_6x6_F_SI1_double(arg0: Annotated[numpy.typing.ArrayLike, numpy.float64])
property T: SmallMatrix_6x6_F_SI1_double
property column_size: int
dot(arg0: SmallMatrix_6x6_F_SI1_double) float[source]
static identity() SmallMatrix_6x6_F_SI1_double[source]
property order: str
prod() float[source]
property row_size: int
set_val(arg0: SupportsFloat | SupportsIndex) SmallMatrix_6x6_F_SI1_double[source]
property size: int
property starting_index: int
sum() float[source]
to_cupy(copy=False, order='F')[source]

Provide a CuPy view into an SmallMatrix.

Note on the order of indices: By default, this is as in AMReX in Fortran contiguous order, indexing as x,y,z. This has performance implications for use in external libraries such as cupy. The order=”C” option will index as z,y,x and perform better with cupy. https://github.com/AMReX-Codes/pyamrex/issues/55#issuecomment-1579610074

Parameters

selfamrex.SmallMatrix_*

A SmallMatrix class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

orderstring, optional

F order (default) or C. C is faster with external libraries.

Returns

cupy.array

A cupy 2-dimensional array.

Raises

ImportError

Raises an exception if cupy is not installed

to_numpy(copy=False, order='F')[source]

Provide a NumPy view into an SmallMatrix.

Note on the order of indices: By default, this is as in AMReX in Fortran contiguous order, indexing as x,y,z. This has performance implications for use in external libraries such as cupy. The order=”C” option will index as z,y,x and perform better with cupy. https://github.com/AMReX-Codes/pyamrex/issues/55#issuecomment-1579610074

Parameters

selfamrex.SmallMatrix_*

A SmallMatrix class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

orderstring, optional

F order (default) or C. C is faster with external libraries.

Returns

np.array

A NumPy 2-dimensional array.

to_xp(copy=False, order='F')[source]

Provide a NumPy or CuPy view into a SmallMatrix, depending on amr.Config.have_gpu .

This function is similar to CuPy’s xp naming suggestion for CPU/GPU agnostic code: https://docs.cupy.dev/en/stable/user_guide/basic.html#how-to-write-cpu-gpu-agnostic-code

Note on the order of indices: By default, this is as in AMReX in Fortran contiguous order, indexing as x,y,z. This has performance implications for use in external libraries such as cupy. The order=”C” option will index as z,y,x and perform better with cupy. https://github.com/AMReX-Codes/pyamrex/issues/55#issuecomment-1579610074

Parameters

selfamrex.SmallMatrix_*

A SmallMatrix class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

orderstring, optional

F order (default) or C. C is faster with external libraries.

Returns

xp.array

A NumPy or CuPy 2-dimensional array.

trace() float[source]
transpose_in_place() SmallMatrix_6x6_F_SI1_double[source]
static zero() SmallMatrix_6x6_F_SI1_double[source]

Utility

class amrex.space3d.ParmParse(prefix: str = '')[source]
add(**kwds)

Helper for @overload to raise when called.

addarr(**kwds)

Helper for @overload to raise when called.

static addfile(arg0: str) None[source]
get_bool(name: str, ival: SupportsInt | SupportsIndex = 0) bool[source]

parses input values

get_int(name: str, ival: SupportsInt | SupportsIndex = 0) int[source]

parses input values

get_real(name: str, ival: SupportsInt | SupportsIndex = 0) float[source]

parses input values

pretty_print_table() None[source]

Write the table in a pretty way to the ostream. If there are duplicates, only the last one is printed.

query_int(name: str, ival: SupportsInt | SupportsIndex = 0) tuple[bool, int][source]

queries input values

remove(arg0: str) int[source]
to_dict() dict[source]

Convert to a nested Python dictionary.

# Example: dump all ParmParse entries to YAML or TOML
import toml
import yaml

pp = amr.ParmParse("").to_dict()
yaml_string = yaml.dump(d)
toml_string = toml.dumps(d)
amrex.space3d.Print(*args, **kwargs)[source]

Wrap amrex::Print() - only the IO processor writes

amrex.space3d.d_decl(x, y, z)[source]

Return a tuple of the three passed elements

amrex.space3d.concatenate(root: str, num: SupportsInt | SupportsIndex, mindigits: SupportsInt | SupportsIndex = 5) str[source]

Builds plotfile name

amrex.space3d.write_single_level_plotfile(plotfilename: str, mf: MultiFab, varnames: Vector_string, geom: Geometry, time: SupportsFloat | SupportsIndex, level_step: SupportsInt | SupportsIndex, versionName: str = 'HyperCLaw-V1.1', levelPrefix: str = 'Level_', mfPrefix: str = 'Cell', extra_dirs: Vector_string = Ellipsis) None[source]

Writes single level plotfile

AmrCore

class amrex.space3d.AmrInfo[source]
blocking_factor(arg0: SupportsInt | SupportsIndex) IntVect3D[source]
check_input: bool
property grid_eff: float
iterate_on_new_grids: bool
max_grid_size(arg0: SupportsInt | SupportsIndex) IntVect3D[source]
property max_level: int
n_error_buf(arg0: SupportsInt | SupportsIndex) IntVect3D[source]
property n_proper: int
ref_ratio(arg0: SupportsInt | SupportsIndex) IntVect3D[source]
refine_grid_layout: bool
refine_grid_layout_dims: IntVect3D
use_fixed_coarse_grids: bool
property use_fixed_upto_level: int
use_new_chop: bool
property verbose: int
class amrex.space3d.AmrMesh[source]
class amrex.space3d.AmrMesh(rb: RealBox, max_level_in: SupportsInt | SupportsIndex, n_cell_in: Vector_int, coord: SupportsInt | SupportsIndex, ref_ratios: Vector_IntVect, is_per: Annotated[collections.abc.Sequence[SupportsInt | SupportsIndex], 'FixedSize(3)'])
property finest_level: int
property max_level: int
ref_ratio(**kwds)

Helper for @overload to raise when called.

property verbose: int

Particles

Particle support is implemented for both legacy (AoS+SoA) and pure SoA particle memory layouts in AMReX. Additional runtime attributes (Real or Int) are always in SoA memory layout.

amrex::StructOfArrays<NReal, NInt, Allocator> is implemented for many numbers of Real and Int arguments, and allocators, e.g.,

amrex::ParticleTile<T_ParticleType, NArrayReal, NArrayInt, Allocator> is implemented for both legacy (AoS+SoA) and pure SoA particle types, many number of Real and Int arguments, and allocators, e.g.,

amrex::ParticleTileData<T_ParticleType, NArrayReal> is implemented for both legacy (AoS+SoA) and pure SoA particle types, many number of Real and Int arguments, e.g.,

amrex::ParticleContainer_impl<ParticleType, T_NArrayReal, T_NArrayInt, Allocator> is implemented for both legacy (AoS+SoA) and pure SoA particle types, many number of Real and Int arguments, and allocators, e.g.,

class amrex.space3d.ParticleContainer_2_1_3_1_default[source]
class amrex.space3d.ParticleContainer_2_1_3_1_default(arg0: Geometry, arg1: DistributionMapping, arg2: BoxArray)
class amrex.space3d.ParticleContainer_2_1_3_1_default(arg0: Vector_Geometry, arg1: Vector_DistributionMapping, arg2: Vector_BoxArray, arg3: Vector_int)
class amrex.space3d.ParticleContainer_2_1_3_1_default(arg0: Vector_Geometry, arg1: Vector_DistributionMapping, arg2: Vector_BoxArray, arg3: Vector_IntVect)
ConstIterator

alias of ParConstIter_2_1_3_1_default

Define(**kwds)

Helper for @overload to raise when called.

Iterator

alias of ParIter_2_1_3_1_default

OK(lev_min: SupportsInt | SupportsIndex = 0, lev_max: SupportsInt | SupportsIndex = -1, nGrow: SupportsInt | SupportsIndex = 0) bool[source]
add_int_comp(**kwds)

Helper for @overload to raise when called.

add_particles(other: ParticleContainer_2_1_3_1_default, local: bool = False) None[source]
add_particles_at_level(particles: ParticleTile_2_1_3_1_default, level: SupportsInt | SupportsIndex, ngrow: SupportsInt | SupportsIndex = 0) None[source]
add_real_comp(**kwds)

Helper for @overload to raise when called.

arena: Arena
property byte_spread: Annotated[list[int], 'FixedSize(3)']
clear_particles() None[source]
const_iterator(*args, level=None)[source]

Create an iterator over all particle tiles

selfamrex.ParticleContainer_*

A ParticleContainer class in pyAMReX

args : deprecated positional argument level : int | str, optional

The MR level. Allowed values are [0:self.finest_level+1) and “all”. If there is more than one MR level, the argument is required.

Iterator over all particle tiles at the specified level.

>>> pc.iterator(level="all")
>>> pc.iterator(level=0)  # only particles on the the coarsest MR level
property finest_level: int
get_int_comp_index(arg0: str) int[source]

Get the Integer SoA index of a component

get_particles(level: SupportsInt | SupportsIndex) dict[tuple[int, int], ParticleTile_2_1_3_1_default][source]
get_real_comp_index(arg0: str) int[source]

Get the ParticleReal SoA index of a component

has_int_comp(arg0: str) bool[source]

Check if a container has an Integer component

has_real_comp(arg0: str) bool[source]

Check if a container has an ParticleReal component

increment(arg0: MultiFab, arg1: SupportsInt | SupportsIndex) None[source]
init_one_per_cell(arg0: SupportsFloat | SupportsIndex, arg1: SupportsFloat | SupportsIndex, arg2: SupportsFloat | SupportsIndex, arg3: ParticleInitType_2_1_3_1) None[source]
init_random(arg0: SupportsInt | SupportsIndex, arg1: SupportsInt | SupportsIndex, arg2: ParticleInitType_2_1_3_1, arg3: bool, arg4: RealBox) None[source]
init_random_per_box(arg0: SupportsInt | SupportsIndex, arg1: SupportsInt | SupportsIndex, arg2: ParticleInitType_2_1_3_1) None[source]
property int_soa_names: list[str]

Get the names for the int SoA components

is_soa_particle: ClassVar[bool] = False
iterator(*args, level=None)[source]

Create an iterator over all particle tiles

selfamrex.ParticleContainer_*

A ParticleContainer class in pyAMReX

args : deprecated positional argument level : int | str, optional

The MR level. Allowed values are [0:self.finest_level+1) and “all”. If there is more than one MR level, the argument is required.

Iterator over all particle tiles at the specified level.

>>> pc.iterator(level="all")
>>> pc.iterator(level=0)  # only particles on the the coarsest MR level
make_alike() ParticleContainer_2_1_3_1_default[source]
num_array_int: ClassVar[int] = 1
num_array_real: ClassVar[int] = 3
property num_int_comps: int

The number of compile-time and runtime int components in SoA

num_local_tiles_at_level(level: SupportsInt | SupportsIndex) int[source]
property num_position_components: int
property num_real_comps: int

The number of compile-time and runtime Real components in SoA

property num_runtime_int_comps: int

The number of runtime Int components in SoA

property num_runtime_real_comps: int

The number of runtime Real components in SoA

num_struct_int: ClassVar[int] = 1
num_struct_real: ClassVar[int] = 2
number_of_particles(only_local: bool = False) int[source]

Return the number of valid particles on all MPI ranks, unless only_local is specified.

number_of_particles_at_level(level: SupportsInt | SupportsIndex, only_valid: bool = True, only_local: bool = False) int[source]
number_of_particles_in_grid(level: SupportsInt | SupportsIndex, only_valid: bool = True, only_local: bool = False) Vector_Long[source]
print_capacity() Annotated[list[int], 'FixedSize(3)'][source]
property real_soa_names: list[str]

Get the names for the Real SoA components

redistribute(lev_min: SupportsInt | SupportsIndex = 0, lev_max: SupportsInt | SupportsIndex = -1, nGrow: SupportsInt | SupportsIndex = 0, local: SupportsInt | SupportsIndex = 0, remove_negative: bool = True) None[source]
remove_particles_at_level(arg0: SupportsInt | SupportsIndex) None[source]
remove_particles_not_at_finestLevel() None[source]
reserve_data() None[source]
resize_data() None[source]
restart(dir: str, file: str) None[source]
restart_checkpoint(dir: str, file: str, is_checkpoint: bool) None[source]
set_soa_compile_time_names(arg0: Sequence[str], arg1: Sequence[str]) None[source]
shrink_t_fit() None[source]
property size: int

Return the number of valid particles on all MPI ranks

sort_particles_by_bin(arg0: IntVect3D) None[source]
sort_particles_by_cell() None[source]
to_df(local=True, comm=None, root_rank=0)[source]

Copy all particles into a pandas.DataFrame

Parameters

selfamrex.ParticleContainer_*

A ParticleContainer class in pyAMReX

localbool

MPI rank-local particles only

commMPI Communicator

if local is False, this defaults to mpi4py.MPI.COMM_WORLD

root_rankMPI root rank to gather to

if local is False, this defaults to 0

Returns

A concatenated pandas.DataFrame with particles from all levels.

Returns None if no particles were found. If local=False, then all ranks but the root_rank will return None.

total_number_of_particles(only_valid: bool = True, only_local: bool = False) int[source]

Return the number of particles (only valid or including invalid) on all MPI ranks, unless only_local is specified.

write_plotfile(dir: str, name: str) None[source]

Likewise for other classes accessible and usable on particle containers:

class amrex.space3d.ParticleInitType_2_1_3_1[source]
property int_array_data: Annotated[list[int], 'FixedSize(1)']
property int_struct_data: Annotated[list[int], 'FixedSize(1)']
is_soa_particle: ClassVar[bool] = False
property real_array_data: Annotated[list[float], 'FixedSize(3)']
property real_struct_data: Annotated[list[float], 'FixedSize(2)']

AoS

This is for the legacy, AoS + SoA particle containers only:

amrex::ArrayOfStructs<T_ParticleType, Allocator> is implemented for many numbers of extra Real and Int arguments, and allocators, e.g.,

class amrex.space3d.ArrayOfStructs_2_1_default[source]
back() Particle_2_1[source]

get back member. Problem!!!!! this is perfo

empty(**kwds)

Helper for @overload to raise when called.

getNumNeighbors() int[source]
numNeighborParticles() int[source]
numParticles() int[source]
numRealParticles() int[source]
numTotalParticles() int[source]
pop_back() None[source]
push_back(arg0: Particle_2_1) None[source]
setNumNeighbors(arg0: SupportsInt | SupportsIndex) None[source]
size() int[source]
static test_sizes() None[source]
to_cupy(copy=False)[source]

Provide CuPy views into a ArrayOfStructs.

Parameters

selfamrex.ArrayOfStructs_*

An ArrayOfStructs class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

Returns

namedtuple

A tuple with real and int components that are each lists of 1D NumPy arrays.

Raises

ImportError

Raises an exception if cupy is not installed

to_host() ArrayOfStructs_2_1_pinned[source]
to_numpy(copy=False)[source]

Provide NumPy views into a ArrayOfStructs.

Parameters

selfamrex.ArrayOfStructs_*

An ArrayOfStructs class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

Returns

namedtuple

A tuple with real and int components that are each lists of 1D NumPy arrays.

to_xp(copy=False)[source]

Provide NumPy or CuPy views into a ArrayOfStructs, depending on amr.Config.have_gpu .

This function is similar to CuPy’s xp naming suggestion for CPU/GPU agnostic code: https://docs.cupy.dev/en/stable/user_guide/basic.html#how-to-write-cpu-gpu-agnostic-code

Parameters

selfamrex.ArrayOfStructs_*

An ArrayOfStructs class in pyAMReX

copybool, optional

Copy the data if true, otherwise create a view (default).

Returns

namedtuple

A tuple with real and int components that are each lists of 1D NumPy or CuPy arrays.

amrex::Particle<T_NReal, T_NInt> is implemented for many numbers of extra Real and Int arguments, e.g.,

class amrex.space3d.Particle_2_1[source]
class amrex.space3d.Particle_2_1(arg0: SupportsFloat | SupportsIndex, arg1: SupportsFloat | SupportsIndex, arg2: SupportsFloat | SupportsIndex)
class amrex.space3d.Particle_2_1(arg0: SupportsFloat | SupportsIndex, arg1: SupportsFloat | SupportsIndex, arg2: SupportsFloat | SupportsIndex, *args)
class amrex.space3d.Particle_2_1(arg0: SupportsFloat | SupportsIndex, arg1: SupportsFloat | SupportsIndex, arg2: SupportsFloat | SupportsIndex, **kwargs)
class amrex.space3d.Particle_2_1(**kwargs)
NInt: ClassVar[int] = 1
NReal: ClassVar[int] = 2
NextID(**kwds)

Helper for @overload to raise when called.

cpu() int[source]
get_idata(**kwds)

Helper for @overload to raise when called.

get_rdata(**kwds)

Helper for @overload to raise when called.

id() int[source]
pos(**kwds)

Helper for @overload to raise when called.

setPos(**kwds)

Helper for @overload to raise when called.

set_idata(**kwds)

Helper for @overload to raise when called.

set_rdata(**kwds)

Helper for @overload to raise when called.

property x: float
property y: float
property z: float

Embedded Boundaries

Embedded boundary (EB) support in pyAMReX is still minimal. To build pyAMReX with EB support, you need to add -DAMReX_EB=ON to CMake build options.

amrex.space3d.EB2_Build(geom: Geometry, required_coarsening_level: SupportsInt | SupportsIndex, max_coarsening_level: SupportsInt | SupportsIndex, ngrow: SupportsInt | SupportsIndex = 4, build_coarse_level_by_coarsening: bool = True, extend_domain_face: bool = True, num_coarsen_opt: SupportsInt | SupportsIndex = 0) None[source]

EB generation

class amrex.space3d.EBFArrayBoxFactory[source]
getVolFrac() MultiFab[source]

Return volume faction MultiFab

amrex.space3d.makeEBFabFactory(geom: Geometry, ba: BoxArray, dm: DistributionMapping, ngrow: Vector_int, support: EBSupport) EBFArrayBoxFactory[source]

Make EBFArrayBoxFactory for given Geometry, BoxArray and DistributionMapping