pyAMReX – Python bindings for AMReX

Version

v0.0.1

Date

Apr 09, 2021

pyAMReX provides python bindings for the AMReX. The primary use-case for this library is to use it within the ExaWind Simulation Environment. The python bindings can be used on Linux or MacOS operating systems, and require Python v3.6 or later.

Installation

Runtime dependencies

Build time dependencies

In addition to the runtime dependencies you will need the following packages installed on your system to build pyAMReX from source

Optionally, you’ll also need the following for development and generating documentation on your own machine

  • pytest for running unit tests

  • Sphinx for generating this documentation

Building from source

To build from source create a new Python environment (either virtualenv or conda env). You can use the requirements.txt in the root directory to install dependencies via pip. Once this step is complete, please execute the following commands to configure and build the python modules

# Clone the git repo
git clone git@github.com:sayerhs/pyamrex.git
cd pyamrex

# Install dependencies
pip install -r requirements.txt

# Build extensions
python setup.py install

Building a development version

If you are developing pyAMReX you should install the package in development mode using the following commands instead of setup.py install

# Build extensions
python setup.py build_ext --inplace

# Install package in develoment mode
pip install -e .

Once the package is installed in editable mode or (development mode), you can execute the build_ext command after editing the Cython files and have the latest version of the compiled extensions available within your environment.

Common build issues

If you are using Anaconda Python (or Conda), please make sure that you install mpi4py via source and that the MPI you used to build these packages are consistent with the ones used to build AMReX. Incompatibilities amongst MPI libraries used to build AMReX can cause spurious memory errors and error messages.

pyAMReX API Reference

AMReX base

class amrex.amrex_base.AMReX

AMReX function wrappers

static finalize()

Call AMReX finalize on exit

static initialize(unicode inp_file=None, list args=None, Comm comm=MPI.COMM_WORLD, unicode exe_name=u'python_driver')

Initialize AMReX library

Parameters
  • inp_file (str) – AMReX input file to be parsed

  • args (list) – A list of “command line arguments” to be passed to ParmParse

  • comm – The MPI communicator instance

  • exe_name (str) – The executable name that AMReX should use as argv[0]

class amrex.amrex_base.Array4Int

AMReX Array4Int interface

This datastructure is usually used from within an MFIter loop

# Assuming BoxArray and DistributionMap are available
ifab = amrex.IMultiFab(ba, dm, 1, 0)
for mfi in ifab:
    bx = mfi.tilebox()
    iarr = ifab.array(mfi)
    iarr[i, j, k, 0] = 10
num_comp

Number of components

size

Size of array

class amrex.amrex_base.Array4Real

AMReX Array4 interface

This datastructure is usually used from within an MFIter loop

# Assuming BoxArray and DistributionMap are available
mfab = amrex.MultiFab(ba, dm, 1, 0)
for mfi in mfab:
    bx = mfi.tilebox()
    marr = mfab.array(mfi)
    marr[i, j, k, 0] = 10
class amrex.amrex_base.Box

Rectangular domain on an Integer Lattice

# Create a Box instance
bx = Box((0, 0, 0), (31, 31, 31))
# Loop over the box
for i, j, k in bx:
    print(i, j, k)
convert(self, CellIndex i, CellIndex j, CellIndex k)

Convert a box to a given type

enclosed_cells(self, int dir=-1)

Convert to a cell centered box

If dir is less than 0 then this method does it for all three directions. If dir is [0, 1, 2] the box is converted only in that direction.

grow(self, int ng)

Grow box by given number of ghost cells

static new(tuple ilo, tuple ihi, tuple itype=None)

Create a new amrex::Box instance

Parameters
  • ilo (tuple) – Three integers defining the low corner

  • ihi (tupe) – Three integers defining the high corner

surrounding_nodes(self, int dir=-1)

Convert to nodal box

If dir is less than 0 then this method does it for all three directions. If dir is [0, 1, 2] the box is converted to nodal only in that direction (i.e., a face box)

hi_vect

Return the high end

is_empty

Flag indicating whether the box is empty

lo_vect

Return the low end

num_pts

Return the number of cells

ok

Flag indicating if the box is ok

class amrex.amrex_base.BoxArray

amrex::BoxArray interface

# Create a new BoxArray instance
bx = Box.new((0,0,0), (127, 127, 127))
ba = BoxArray.new(ba)

# Chop up the boxes by setting max grid size
ba.max_size(64)
clear(self)

Clear the box array

enclosed_cells(self, int dir=-1)

Return a cell-centered box array

grow(self, int n)

Grow all boxes by given number of ghost cells

max_size(self, int block_size)

Set the max grid size and chop up the BoxArray

static new(Box bx)
resize(self, Long blen)
surrounding_nodes(self, int dir=-1)

Return a nodal box array

empty

Flag indicating if the BoxArray is empty

ok

Flag indicating if the BoxArray is initialized properly

size

Number of boxes in this array

class amrex.amrex_base.DistributionMapping

amrex::DistributionMapping instance

# Create a new BoxArray instance
bx = Box.new((0,0,0), (127, 127, 127))
ba = BoxArray.new(ba)

# Chop up the boxes by setting max grid size
ba.max_size(64)

# Create a distribution map based on the boxarray
dm = DistributionMapping.new(ba)
static new(BoxArray ba, int nprocs=0)
class amrex.amrex_base.Geometry

Rectangular problem domain geometry

is_periodic(self, int dir)

Flag indicating whether domain is periodic in given direction

period(self, int dir)

Periodic in given direction

cell_size

Return the cell size

inv_cell_size

Inverse cell size

is_all_periodic

Flag indicating if all sides are periodic

prob_domain

Return the RealBox for problem domain

prob_hi

Return the high end of the problem domain

prob_lo

Return the low end of the problem domain

prob_size

Return the size of the domain

class amrex.amrex_base.IMultiFab

AMReX IMultiFab wrapper

add(self, IMultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)
array(self, MFIter mfi)
copy(self, IMultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)
divide(self, IMultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)
invert(self, int val, int comp, int ncomp)
is_cell_centered(self)
is_nodal(self, int dir=-1)
max(self, int ncomp=0, int nghost=0)

Return the maximum value for a given component

min(self, int ncomp=0, int nghost=0)

Return the minimum value for a given component

mult(self, int val, int comp, int ncomp)
multiply(self, IMultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)
static new(BoxArray ba, DistributionMapping dm, int ncomp, int nghost)
plus(self, int val, int comp, int ncomp)
set_boundary(self, int val)
set_val(self, int val, int start_comp, int num_comp, int num_ghost=0)

Initialize field to a given value

subtract(self, IMultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)
class amrex.amrex_base.MFIter

MFIter wrapper

fabbox(self)
growntilebox(self)
tilebox(self)
validbox(self)
class amrex.amrex_base.MultiFab

amrex::MultiFab interface

# Create a new BoxArray instance
bx = Box.new((0,0,0), (127, 127, 127))
ba = BoxArray.new(ba)

# Chop up the boxes by setting max grid size
ba.max_size(64)

# Create a distribution map based on the boxarray
dm = DistributionMapping.new(ba)

# Create a MultiFab instance
num_components = 3
num_ghost = 0
mfab = MultiFab.new(ba, dm, num_components, num_ghost)
# Initialize the field to [10.0, 20.0, 30.0]
mfab.set_val(10.0, 0, 1)
mfab.set_val(20.0, 1, 1)
mfab.set_val(30.0, 2, 1)

# Perform some operations
mfab.plus(10.0, 0, 1)
mfab.mult(2.0, 1, 1)
mfab.plus(5.0, 2, 1)

# Looping
for mfi in mfab:
    bx = mfi.tilebox()
    marr = mfab.array(mfi)

    for i, j, k in bx:
        mar[i, j, k, 0] = 10.0 * i
        mar[i, j, k, 1] = 10.0 * j
        mar[i, j, k, 2] = 10.0 * k
abs(self, int comp, int ncomp)

x[:, comp:comp+ncomp] = abs(x[:, comp:comp+ncomp])

add(self, MultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)

dst[:, dcomp:dcomp+ncomp] += src[:, scomp:scomp+ncomp]

array(self, MFIter mfi)
copy(self, MultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)

dst[:, dcomp:dcomp+ncomp] = src[:, scomp:scomp+ncomp]

divide(self, MultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)

dst[:, dcomp:dcomp+ncomp] /= src[:, scomp:scomp+ncomp]

invert(self, Real val, int comp, int ncomp)

x[:, comp:comp+ncomp] = (val / x[:, comp:comp+ncomp])

is_cell_centered(self)
is_nodal(self, int dir=-1)
lin_comb(self, Real aval, MultiFab x, int xcomp, Real bval, MultiFab y, int ycomp, int dcomp, int ncomp, int nghost=0)

dst = a * x + b * y

max(self, int ncomp=0, int nghost=0)

Return the maximum value for a given component

min(self, int ncomp=0, int nghost=0)

Return the minimum value for a given component

mult(self, Real val, int comp, int ncomp)

x[:, comp:comp+ncomp] *= val

multiply(self, MultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)

dst[:, dcomp:dcomp+ncomp] *= src[:, scomp:scomp+ncomp]

static new(BoxArray ba, DistributionMapping dm, int ncomp, int nghost)
plus(self, Real val, int comp, int ncomp)

x[:, comp:comp+ncomp] += val

saxpy(self, Real aval, MultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)

dst[:, dcomp:dcomp+ncomp] += aval * src[:, scomp:scomp+ncomp]

set_boundary(self, Real val)
set_val(self, Real val, int start_comp, int num_comp, int num_ghost=0)

Initialize field to a given value

subtract(self, MultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)

dst[:, dcomp:dcomp+ncomp] -= src[:, scomp:scomp+ncomp]

swap(self, MultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)
xpay(self, Real aval, MultiFab src, int scomp=0, int dcomp=0, int ncomp=1, int nghost=0)

dst[:, dcomp:dcomp+ncomp] = src[:, scomp:scomp+ncomp] + aval * dst[:, dcomp:dcomp+ncomp]

class amrex.amrex_base.ParmParse

AMReX parameter parser

Example usage:

pp = amrex.ParmParse("ABL")
src_terms = pp.query("source_terms")
tvals = pp.query_real("temperature_values")
add(self, unicode key, val)

Add a single value to the database

Parameters
  • key (str) – Keyword to access the value

  • val – An int, real, or string value

add_arr(self, unicode key, list values)

Add a list of values corresponding to a lookup key in the table

Parameters
  • key (str) – Keyword to access the value

  • valuse – A list of int, real, or string values

get(self, unicode key)

Return value as string

get_impl(self, unicode key, func)
Parameters
  • key (str) – Keyword to lookup

  • func – Conversion function

get_int(self, unicode key)

Return value as string

get_real(self, unicode key)

Return value as string

query(self, unicode key)

Return if a value exists in the input file else None

Parameters

key (str) – Keyword to lookup

Returns

str or None

query_int(self, unicode key)

Return if an integer value exists or else return None

Parameters

key (str) – Keyword to lookup

Returns

int or None

query_real(self, unicode key)

Return if a value exists in the input file else None

Parameters

key (str) – Keyword to lookup

Returns

real or None

class amrex.amrex_base.RealBox

Domain box with dimensions

length(self, int dir)

Return length in given direction

static new(tuple xlo, tuple xhi)

Create a new C++ amrex::RealBox instance

Parameters
  • xlo (tuple) – Coordinates of the lower corner

  • xhi (tuple) – Coordinates of the upper corner

hi

Return the upper corner coordinats as a NumPy array

lo

Return lower corner coordinates as a NumPy array

AMReX AmrCore bindings

class amrex.amrex_core.AmrCore

AmrCore wrapper

boxarray(self, lev=None)

Return BoxArray

dmap(self, lev=None)

Return distribution mapping

geom(self, lev=None)

Return geometry instance

finest_level

Current finest level

max_level

Return the maximum level possible

num_levels

Total number of levels

verbose

Return verbosity

Indices and tables