Developer Guide

A Brief Overview of the PyFR Framework

Where to Start

The symbolic link pyfr.scripts.pyfr points to the script pyfr.scripts.main, which is where it all starts! Specifically, the function process_run calls the function _process_common, which in turn calls the function get_solver, returning an Integrator – a composite of a Controller and a Stepper. The Integrator has a method named run, which is then called to run the simulation.

Controller

A Controller acts to advance the simulation in time. Specifically, a Controller has a method named advance_to which advances a System to a specified time. There are three types of Controller available in PyFR 1.5.0:

class pyfr.integrators.std.controllers.StdNoneController(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_accept_step(dt, idxcurr, err=None)
_add(*args)
_controller_needs_errest
_get_axnpby_kerns(n, subdims=None)
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_reject_step(dt, idxold, err=None)
_stepper_has_errest
_stepper_nfevals
_stepper_nregs
_stepper_order
advance_to(t)[source]
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
controller_name = 'none'
formulation = 'std'
nsteps
run()
soln
step(t, dt)
class pyfr.integrators.std.controllers.StdPIController(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_accept_step(dt, idxcurr, err=None)
_add(*args)
_controller_needs_errest
_errest(x, y, z)[source]
_get_axnpby_kerns(n, subdims=None)
_get_errest_kerns()[source]
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_reject_step(dt, idxold, err=None)
_stepper_has_errest
_stepper_nfevals
_stepper_nregs
_stepper_order
advance_to(t)[source]
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
controller_name = 'pi'
formulation = 'std'
nsteps
run()
soln
step(t, dt)
class pyfr.integrators.dual.controllers.DualNoneController(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_accept_step(dt, idxcurr)
_add(*args)
_dual_time_source()
_get_axnpby_kerns(n, subdims=None)
_get_errest_kerns()[source]
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_resid(x, y)[source]
_source_regidx
_stepper_nfevals
_stepper_nregs
_stepper_order
_stepper_regidx
advance_to(t)[source]
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
controller_name = 'none'
finalise_step(currsoln)
formulation = 'dual'
nsteps
run()
soln
step(t, dt)

Types of Controller are related via the following inheritance diagram:

Inheritance diagram of pyfr.integrators.std.controllers, pyfr.integrators.dual.controllers

Stepper

A Stepper acts to advance the simulation by a single time-step. Specifically, a Stepper has a method named step which advances a System by a single time-step. There are 11 types of Stepper available in PyFR 1.5.0:

class pyfr.integrators.std.steppers.StdEulerStepper(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_add(*args)
_controller_needs_errest
_get_axnpby_kerns(n, subdims=None)
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_stepper_has_errest
_stepper_nfevals
_stepper_nregs
_stepper_order
advance_to(t)
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
formulation = 'std'
nsteps
run()
soln
step(t, dt)[source]
stepper_name = 'euler'
class pyfr.integrators.std.steppers.StdRK4Stepper(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_add(*args)
_controller_needs_errest
_get_axnpby_kerns(n, subdims=None)
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_stepper_has_errest
_stepper_nfevals
_stepper_nregs
_stepper_order
advance_to(t)
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
formulation = 'std'
nsteps
run()
soln
step(t, dt)[source]
stepper_name = 'rk4'
class pyfr.integrators.std.steppers.StdRK34Stepper(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_add(*args)
_controller_needs_errest
_get_axnpby_kerns(n, subdims=None)
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_stepper_has_errest
_stepper_nfevals
_stepper_nregs
_stepper_order
a = [0.32416573882874605, 0.5570978645055429, -0.08605491431272755]
advance_to(t)
b = [0.10407986927510238, 0.6019391368822611, 2.9750900268840206, -2.681109033041384]
bhat = [0.3406814840808433, 0.09091523008632837, 2.866496742725443, -2.298093456892615]
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
formulation = 'std'
nsteps
run()
soln
step(t, dt)
stepper_name = 'rk34'
class pyfr.integrators.std.steppers.StdRK45Stepper(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_add(*args)
_controller_needs_errest
_get_axnpby_kerns(n, subdims=None)
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_stepper_has_errest
_stepper_nfevals
_stepper_nregs
_stepper_order
a = [0.22502245872571303, 0.5440433129514047, 0.14456824349399464, 0.7866643421983568]
advance_to(t)
b = [0.05122930664033915, 0.3809548257264019, -0.3733525963923833, 0.5925012850263623, 0.34866717899927996]
bhat = [0.13721732210321927, 0.19188076232938728, -0.2292067211595315, 0.6242946765438954, 0.27581396018302956]
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
formulation = 'std'
nsteps
run()
soln
step(t, dt)
stepper_name = 'rk45'
class pyfr.integrators.std.steppers.StdTVDRK3Stepper(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_add(*args)
_controller_needs_errest
_get_axnpby_kerns(n, subdims=None)
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_stepper_has_errest
_stepper_nfevals
_stepper_nregs
_stepper_order
advance_to(t)
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
formulation = 'std'
nsteps
run()
soln
step(t, dt)[source]
stepper_name = 'tvd-rk3'
class pyfr.integrators.dual.steppers.DualBDF2Stepper(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_add(*args)
_dual_time_source
_get_axnpby_kerns(n, subdims=None)
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_source_regidx
_stepper_nfevals
_stepper_nregs
_stepper_order
_stepper_regidx
advance_to(t)
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
finalise_step(currsoln)
formulation = 'dual'
nsteps
run()
soln
step(t, dt)
stepper_name = 'bdf2'
class pyfr.integrators.dual.steppers.DualBDF3Stepper(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_add(*args)
_dual_time_source
_get_axnpby_kerns(n, subdims=None)
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_source_regidx
_stepper_nfevals
_stepper_nregs
_stepper_order
_stepper_regidx
advance_to(t)
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
finalise_step(currsoln)
formulation = 'dual'
nsteps
run()
soln
step(t, dt)
stepper_name = 'bdf3'
class pyfr.integrators.dual.steppers.DualBackwardEulerStepper(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_add(*args)
_dual_time_source
_get_axnpby_kerns(n, subdims=None)
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_source_regidx
_stepper_nfevals
_stepper_nregs
_stepper_order
_stepper_regidx
advance_to(t)
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
finalise_step(currsoln)
formulation = 'dual'
nsteps
run()
soln
step(t, dt)
stepper_name = 'backward-euler'
class pyfr.integrators.dual.pseudosteppers.DualPseudoRK4Stepper(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_add(*args)
_add_with_dts(*args, c)
_dual_time_source()
_get_axnpby_kerns(n, subdims=None)
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_pseudo_stepper_nregs
_pseudo_stepper_order
_source_regidx
_stepper_nfevals
_stepper_nregs
_stepper_order
_stepper_regidx
advance_to(t)
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
finalise_step(currsoln)
formulation = 'dual'
nsteps
pseudo_stepper_name = 'rk4'
run()
soln
step(t, dt, dtau)[source]
class pyfr.integrators.dual.pseudosteppers.DualPseudoTVDRK3Stepper(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_add(*args)
_add_with_dts(*args, c)
_dual_time_source()
_get_axnpby_kerns(n, subdims=None)
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_pseudo_stepper_nregs
_pseudo_stepper_order
_source_regidx
_stepper_nfevals
_stepper_nregs
_stepper_order
_stepper_regidx
advance_to(t)
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
finalise_step(currsoln)
formulation = 'dual'
nsteps
pseudo_stepper_name = 'tvd-rk3'
run()
soln
step(t, dt, dtau)[source]
class pyfr.integrators.dual.pseudosteppers.DualPseudoEulerStepper(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_add(*args)
_add_with_dts(*args, c)
_dual_time_source()
_get_axnpby_kerns(n, subdims=None)
_get_gndofs()
_get_kernels(name, nargs, **kwargs)
_get_plugins()
_get_reg_banks(nreg)
_prepare_reg_banks(*bidxes)
_pseudo_stepper_nregs
_pseudo_stepper_order
_source_regidx
_stepper_nfevals
_stepper_nregs
_stepper_order
_stepper_regidx
advance_to(t)
call_plugin_dt(dt)
cfgmeta
collect_stats(stats)
finalise_step(currsoln)
formulation = 'dual'
nsteps
pseudo_stepper_name = 'euler'
run()
soln
step(t, dt, dtau)[source]

Types of Stepper are related via the following inheritance diagram:

Inheritance diagram of pyfr.integrators.dual.steppers, pyfr.integrators.dual.pseudosteppers, pyfr.integrators.std.steppers

System

A System holds information/data for the system, including Elements, Interfaces, and the Backend with which the simulation is to run. A System has a method named rhs, which obtains the divergence of the flux (the ‘right-hand-side’) at each solution point. The method rhs invokes various kernels which have been pre-generated and loaded into queues. A System also has a method named _gen_kernels which acts to generate all the kernels required by a particular System. A kernel is an instance of a ‘one-off’ class with a method named run that implements the required kernel functionality. Individual kernels are produced by a kernel provider. PyFR 1.5.0 has various types of kernel provider. A Pointwise Kernel Provider produces point-wise kernels such as Riemann solvers and flux functions etc. These point-wise kernels are specified using an in-built platform-independent templating language derived from Mako, henceforth referred to as PyFR-Mako. There are two types of System available in PyFR 1.5.0:

class pyfr.solvers.euler.system.EulerSystem(backend, rallocs, mesh, initsoln, nreg, cfg)[source]
_abc_impl = <_abc_data object>
_gen_kernels(eles, iint, mpiint, bcint)
_gen_queues()
_load_bc_inters(rallocs, mesh, elemap)
_load_eles(rallocs, mesh, initsoln, nreg, nonce)
_load_int_inters(rallocs, mesh, elemap)
_load_mpi_inters(rallocs, mesh, elemap)
_nonce_seq = count(0)
_nqueues = 2
bbcinterscls

alias of pyfr.solvers.euler.inters.EulerBaseBCInters

ele_scal_upts(idx)
elementscls

alias of pyfr.solvers.euler.elements.EulerElements

filt(uinoutbank)
intinterscls

alias of pyfr.solvers.euler.inters.EulerIntInters

mpiinterscls

alias of pyfr.solvers.euler.inters.EulerMPIInters

name = 'euler'
rhs(t, uinbank, foutbank)
class pyfr.solvers.navstokes.system.NavierStokesSystem(backend, rallocs, mesh, initsoln, nreg, cfg)[source]
_abc_impl = <_abc_data object>
_gen_kernels(eles, iint, mpiint, bcint)
_gen_queues()
_load_bc_inters(rallocs, mesh, elemap)
_load_eles(rallocs, mesh, initsoln, nreg, nonce)
_load_int_inters(rallocs, mesh, elemap)
_load_mpi_inters(rallocs, mesh, elemap)
_nonce_seq = count(0)
_nqueues = 2
bbcinterscls

alias of pyfr.solvers.navstokes.inters.NavierStokesBaseBCInters

ele_scal_upts(idx)
elementscls

alias of pyfr.solvers.navstokes.elements.NavierStokesElements

filt(uinoutbank)
intinterscls

alias of pyfr.solvers.navstokes.inters.NavierStokesIntInters

mpiinterscls

alias of pyfr.solvers.navstokes.inters.NavierStokesMPIInters

name = 'navier-stokes'
rhs(t, uinbank, foutbank)

Types of System are related via the following inheritance diagram:

Inheritance diagram of pyfr.solvers.navstokes.system, pyfr.solvers.euler.system

Elements

An Elements holds information/data for a group of elements. There are two types of Elements available in PyFR 1.5.0:

class pyfr.solvers.euler.elements.EulerElements(basiscls, eles, cfg)[source]
_abc_impl = <_abc_data object>
_gen_pnorm_fpts()
_mag_pnorm_fpts = None
_norm_pnorm_fpts = None
_ploc_in_src_exprs = None
_scratch_bufs
_smats_djacs_mpts = None
_soln_in_src_exprs = None
_src_exprs = None
_srtd_face_fpts = None
static con_to_pri(cons, cfg)
convarmap = {2: ['rho', 'rhou', 'rhov', 'E'], 3: ['rho', 'rhou', 'rhov', 'rhow', 'E']}
dualcoeffs = {2: ['rho', 'rhou', 'rhov', 'E'], 3: ['rho', 'rhou', 'rhov', 'rhow', 'E']}
formulations = ['std', 'dual']
get_mag_pnorms(eidx, fidx)
get_mag_pnorms_for_inter(eidx, fidx)
get_norm_pnorms(eidx, fidx)
get_norm_pnorms_for_inter(eidx, fidx)
get_ploc_for_inter(eidx, fidx)
get_scal_fpts_for_inter(eidx, fidx)
get_vect_fpts_for_inter(eidx, fidx)
opmat(expr)
ploc_at(name)
ploc_at_np(name)
plocfpts = None
static pri_to_con(pris, cfg)
privarmap = {2: ['rho', 'u', 'v', 'p'], 3: ['rho', 'u', 'v', 'w', 'p']}
rcpdjac_at(name)
rcpdjac_at_np(name)
set_backend(backend, nscalupts, nonce)[source]
set_ics_from_cfg()
set_ics_from_soln(solnmat, solncfg)
smat_at(name)
smat_at_np(name)
visvarmap = {2: {'density': ['rho'], 'pressure': ['p'], 'velocity': ['u', 'v']}, 3: {'density': ['rho'], 'pressure': ['p'], 'velocity': ['u', 'v', 'w']}}
class pyfr.solvers.navstokes.elements.NavierStokesElements(basiscls, eles, cfg)[source]
_abc_impl = <_abc_data object>
_gen_pnorm_fpts()
_mag_pnorm_fpts = None
_norm_pnorm_fpts = None
_ploc_in_src_exprs = None
_scratch_bufs
_smats_djacs_mpts = None
_soln_in_src_exprs = None
_src_exprs = None
_srtd_face_fpts = None
static con_to_pri(cons, cfg)
convarmap = {2: ['rho', 'rhou', 'rhov', 'E'], 3: ['rho', 'rhou', 'rhov', 'rhow', 'E']}
dualcoeffs = {2: ['rho', 'rhou', 'rhov', 'E'], 3: ['rho', 'rhou', 'rhov', 'rhow', 'E']}
formulations = ['std', 'dual']
get_artvisc_fpts_for_inter(eidx, fidx)
get_mag_pnorms(eidx, fidx)
get_mag_pnorms_for_inter(eidx, fidx)
get_norm_pnorms(eidx, fidx)
get_norm_pnorms_for_inter(eidx, fidx)
get_ploc_for_inter(eidx, fidx)
get_scal_fpts_for_inter(eidx, fidx)
get_vect_fpts_for_inter(eidx, fidx)
opmat(expr)
ploc_at(name)
ploc_at_np(name)
plocfpts = None
static pri_to_con(pris, cfg)
privarmap = {2: ['rho', 'u', 'v', 'p'], 3: ['rho', 'u', 'v', 'w', 'p']}
rcpdjac_at(name)
rcpdjac_at_np(name)
set_backend(backend, nscalupts, nonce)[source]
set_ics_from_cfg()
set_ics_from_soln(solnmat, solncfg)
shockvar = 'rho'
smat_at(name)
smat_at_np(name)
visvarmap = {2: {'density': ['rho'], 'pressure': ['p'], 'velocity': ['u', 'v']}, 3: {'density': ['rho'], 'pressure': ['p'], 'velocity': ['u', 'v', 'w']}}

Types of Elements are related via the following inheritance diagram:

Inheritance diagram of pyfr.solvers.navstokes.elements, pyfr.solvers.euler.elements

Interfaces

An Interfaces holds information/data for a group of interfaces. There are four types of (non-boundary) Interfaces available in PyFR 1.5.0:

class pyfr.solvers.euler.inters.EulerIntInters(*args, **kwargs)[source]
_const_mat(inter, meth)
_gen_perm(lhs, rhs)
_scal_view(inter, meth)
_scal_xchg_view(inter, meth)
_vect_view(inter, meth)
_vect_xchg_view(inter, meth)
_view(inter, meth, vshape=())
_xchg_view(inter, meth, vshape=())
class pyfr.solvers.euler.inters.EulerMPIInters(*args, **kwargs)[source]
MPI_TAG = 2314
_const_mat(inter, meth)
_scal_view(inter, meth)
_scal_xchg_view(inter, meth)
_vect_view(inter, meth)
_vect_xchg_view(inter, meth)
_view(inter, meth, vshape=())
_xchg_view(inter, meth, vshape=())
class pyfr.solvers.navstokes.inters.NavierStokesIntInters(be, lhs, rhs, elemap, cfg)[source]
_const_mat(inter, meth)
_gen_perm(lhs, rhs)
_scal_view(inter, meth)
_scal_xchg_view(inter, meth)
_vect_view(inter, meth)
_vect_xchg_view(inter, meth)
_view(inter, meth, vshape=())
_xchg_view(inter, meth, vshape=())
class pyfr.solvers.navstokes.inters.NavierStokesMPIInters(be, lhs, rhsrank, rallocs, elemap, cfg)[source]
MPI_TAG = 2314
_const_mat(inter, meth)
_scal_view(inter, meth)
_scal_xchg_view(inter, meth)
_vect_view(inter, meth)
_vect_xchg_view(inter, meth)
_view(inter, meth, vshape=())
_xchg_view(inter, meth, vshape=())

Types of (non-boundary) Interfaces are related via the following inheritance diagram:

Inheritance diagram of pyfr.solvers.navstokes.inters.NavierStokesMPIInters, pyfr.solvers.navstokes.inters.NavierStokesIntInters, pyfr.solvers.euler.inters.EulerMPIInters, pyfr.solvers.euler.inters.EulerIntInters

Backend

A Backend holds information/data for a backend. There are four types of Backend available in PyFR 1.5.0:

class pyfr.backends.cuda.base.CUDABackend(cfg)[source]
_abc_impl = <_abc_data object>
_malloc_impl(nbytes)[source]
alias(obj, aobj)
commit()
const_matrix(initval, extent=None, tags={})
kernel(name, *args, **kwargs)
lookup = None
malloc(obj, extent)
matrix(ioshape, initval=None, extent=None, aliases=None, tags={})
matrix_bank(mats, initbank=0, tags={})
matrix_rslice(mat, p, q)
name = 'cuda'
queue()
runall(sequence)
view(matmap, rmap, cmap, rstridemap=None, vshape=(), tags={})
xchg_matrix(ioshape, initval=None, extent=None, aliases=None, tags={})
xchg_matrix_for_view(view, tags={})
xchg_view(matmap, rmap, cmap, rstridemap=None, vshape=(), tags={})
class pyfr.backends.mic.base.MICBackend(cfg)[source]
_abc_impl = <_abc_data object>
_malloc_impl(nbytes)[source]
alias(obj, aobj)
commit()
const_matrix(initval, extent=None, tags={})
kernel(name, *args, **kwargs)
lookup = None
malloc(obj, extent)
matrix(ioshape, initval=None, extent=None, aliases=None, tags={})
matrix_bank(mats, initbank=0, tags={})
matrix_rslice(mat, p, q)
name = 'mic'
queue()
runall(sequence)
view(matmap, rmap, cmap, rstridemap=None, vshape=(), tags={})
xchg_matrix(ioshape, initval=None, extent=None, aliases=None, tags={})
xchg_matrix_for_view(view, tags={})
xchg_view(matmap, rmap, cmap, rstridemap=None, vshape=(), tags={})
class pyfr.backends.opencl.base.OpenCLBackend(cfg)[source]
_abc_impl = <_abc_data object>
_malloc_impl(nbytes)[source]
alias(obj, aobj)
commit()
const_matrix(initval, extent=None, tags={})
kernel(name, *args, **kwargs)
lookup = None
malloc(obj, extent)
matrix(ioshape, initval=None, extent=None, aliases=None, tags={})
matrix_bank(mats, initbank=0, tags={})
matrix_rslice(mat, p, q)
name = 'opencl'
queue()
runall(sequence)
view(matmap, rmap, cmap, rstridemap=None, vshape=(), tags={})
xchg_matrix(ioshape, initval=None, extent=None, aliases=None, tags={})
xchg_matrix_for_view(view, tags={})
xchg_view(matmap, rmap, cmap, rstridemap=None, vshape=(), tags={})
class pyfr.backends.openmp.base.OpenMPBackend(cfg)[source]
_abc_impl = <_abc_data object>
_malloc_impl(nbytes)[source]
alias(obj, aobj)
commit()
const_matrix(initval, extent=None, tags={})
kernel(name, *args, **kwargs)
lookup = None
malloc(obj, extent)
matrix(ioshape, initval=None, extent=None, aliases=None, tags={})
matrix_bank(mats, initbank=0, tags={})
matrix_rslice(mat, p, q)
name = 'openmp'
queue()
runall(sequence)
view(matmap, rmap, cmap, rstridemap=None, vshape=(), tags={})
xchg_matrix(ioshape, initval=None, extent=None, aliases=None, tags={})
xchg_matrix_for_view(view, tags={})
xchg_view(matmap, rmap, cmap, rstridemap=None, vshape=(), tags={})

Types of Backend are related via the following inheritance diagram:

Inheritance diagram of pyfr.backends.cuda.base, pyfr.backends.opencl.base, pyfr.backends.openmp.base, pyfr.backends.mic.base

Pointwise Kernel Provider

A Pointwise Kernel Provider produces point-wise kernels. Specifically, a Pointwise Kernel Provider has a method named register, which adds a new method to an instance of a Pointwise Kernel Provider. This new method, when called, returns a kernel. A kernel is an instance of a ‘one-off’ class with a method named run that implements the required kernel functionality. The kernel functionality itself is specified using PyFR-Mako. Hence, a Pointwise Kernel Provider also has a method named _render_kernel, which renders PyFR-Mako into low-level platform-specific code. The _render_kernel method first sets the context for Mako (i.e. details about the Backend etc.) and then uses Mako to begin rendering the PyFR-Mako specification. When Mako encounters a pyfr:kernel an instance of a Kernel Generator is created, which is used to render the body of the pyfr:kernel. There are four types of Pointwise Kernel Provider available in PyFR 1.5.0:

class pyfr.backends.cuda.provider.CUDAPointwiseKernelProvider(backend)[source]
_abc_impl = <_abc_data object>
_build_arglst(dims, argn, argt, argdict)
_build_kernel(name, src, argtypes)
_instantiate_kernel(dims, fun, arglst)[source]
_render_kernel(name, mod, tplargs)
kernel_generator_cls

alias of pyfr.backends.cuda.generator.CUDAKernelGenerator

register(mod)
class pyfr.backends.mic.provider.MICPointwiseKernelProvider(backend)[source]
_abc_impl = <_abc_data object>
_build_arglst(dims, argn, argt, argdict)
_build_kernel(name, src, argtypes, restype=None)
_instantiate_kernel(dims, fun, arglst)[source]
_render_kernel(name, mod, tplargs)
kernel_generator_cls

alias of pyfr.backends.mic.generator.MICKernelGenerator

register(mod)
class pyfr.backends.opencl.provider.OpenCLPointwiseKernelProvider(backend)[source]
_abc_impl = <_abc_data object>
_build_arglst(dims, argn, argt, argdict)
_build_kernel(name, src, argtypes)
_instantiate_kernel(dims, fun, arglst)[source]
_render_kernel(name, mod, tplargs)
kernel_generator_cls

alias of pyfr.backends.opencl.generator.OpenCLKernelGenerator

register(mod)
class pyfr.backends.openmp.provider.OpenMPPointwiseKernelProvider(backend)[source]
_abc_impl = <_abc_data object>
_build_arglst(dims, argn, argt, argdict)
_build_kernel(name, src, argtypes, restype=None)
_instantiate_kernel(dims, fun, arglst)[source]
_render_kernel(name, mod, tplargs)
kernel_generator_cls

alias of pyfr.backends.openmp.generator.OpenMPKernelGenerator

register(mod)

Types of Pointwise Kernel Provider are related via the following inheritance diagram:

Inheritance diagram of pyfr.backends.openmp.provider, pyfr.backends.cuda.provider, pyfr.backends.opencl.provider, pyfr.backends.mic.provider, pyfr.backends.base.kernels.BasePointwiseKernelProvider

Kernel Generator

A Kernel Generator renders the PyFR-Mako in a pyfr:kernel into low-level platform-specific code. Specifically, a Kernel Generator has a method named render, which applies Backend specific regex and adds Backend specific ‘boiler plate’ code to produce the low-level platform-specific source – which is compiled, linked, and loaded. There are four types of Kernel Generator available in PyFR 1.5.0:

class pyfr.backends.cuda.generator.CUDAKernelGenerator(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_deref_arg_array_1d(arg)
_deref_arg_array_2d(arg)
_deref_arg_view(arg)
_render_body(body)
_render_spec()[source]
argspec()
needs_ldim(arg)
render()[source]
class pyfr.backends.mic.generator.MICKernelGenerator(name, ndim, args, body, fpdtype)[source]
_abc_impl = <_abc_data object>
_deref_arg_array_1d(arg)
_deref_arg_array_2d(arg)
_deref_arg_view(arg)
_render_body(body)
_render_spec_unpack()[source]
argspec()
needs_ldim(arg)
render()[source]
class pyfr.backends.opencl.generator.OpenCLKernelGenerator(*args, **kwargs)[source]
_abc_impl = <_abc_data object>
_deref_arg_array_1d(arg)
_deref_arg_array_2d(arg)
_deref_arg_view(arg)
_render_body(body)
_render_spec()[source]
argspec()
needs_ldim(arg)
render()[source]
class pyfr.backends.openmp.generator.OpenMPKernelGenerator(name, ndim, args, body, fpdtype)[source]
_abc_impl = <_abc_data object>
_deref_arg_array_1d(arg)
_deref_arg_array_2d(arg)
_deref_arg_view(arg)
_render_body(body)
_render_spec()[source]
argspec()
needs_ldim(arg)
render()[source]

Types of Kernel Generator are related via the following inheritance diagram:

Inheritance diagram of pyfr.backends.cuda.generator.CUDAKernelGenerator, pyfr.backends.opencl.generator.OpenCLKernelGenerator, pyfr.backends.openmp.generator.OpenMPKernelGenerator, pyfr.backends.mic.generator.MICKernelGenerator

PyFR-Mako

PyFR-Mako Kernels

PyFR-Mako kernels are specifications of point-wise functionality that can be invoked directly from within PyFR. They are opened with a header of the form:

<%pyfr:kernel name='kernel-name' ndim='data-dimensionality' [argument-name='argument-intent argument-attribute argument-data-type' ...]>

where

  1. kernel-name — name of kernel

    string

  2. data-dimensionality — dimensionality of data

    int

  3. argument-name — name of argument

    string

  4. argument-intent — intent of argument

    in | out | inout

  5. argument-attribute — attribute of argument

    mpi | scalar | view

  6. argument-data-type — data type of argument

    string

and are closed with a footer of the form:

</%pyfr:kernel>

PyFR-Mako Macros

PyFR-Mako macros are specifications of point-wise functionality that cannot be invoked directly from within PyFR, but can be embedded into PyFR-Mako kernels. PyFR-Mako macros can be viewed as building blocks for PyFR-mako kernels. They are opened with a header of the form:

<%pyfr:macro name='macro-name' params='[parameter-name, ...]'>

where

  1. macro-name — name of macro

    string

  2. parameter-name — name of parameter

    string

and are closed with a footer of the form:

</%pyfr:macro>

PyFR-Mako macros are embedded within a kernel using an expression of the following form:

${pyfr.expand('macro-name', ['parameter-name', ...])};

where

  1. macro-name — name of the macro

    string

  2. parameter-name — name of parameter

    string

Syntax

Basic Functionality

Basic functionality can be expressed using a restricted subset of the C programming language. Specifically, use of the following is allowed:

  1. +,-,*,/ — basic arithmetic
  2. sin, cos, tan — basic trigonometric functions
  3. exp — exponential
  4. pow — power
  5. fabs — absolute value
  6. output = ( condition ? satisfied : unsatisfied ) — ternary if
  7. min — minimum
  8. max — maximum

However, conditional if statements, as well as for/while loops, are not allowed.

Expression Substitution

Mako expression substitution can be used to facilitate PyFR-Mako kernel specification. A Python expression expression prescribed thus ${expression} is substituted for the result when the PyFR-Mako kernel specification is interpreted at runtime.

Example:

E = s[${ndims - 1}]

Conditionals

Mako conditionals can be used to facilitate PyFR-Mako kernel specification. Conditionals are opened with % if condition: and closed with % endif. Note that such conditionals are evaluated when the PyFR-Mako kernel specification is interpreted at runtime, they are not embedded into the low-level kernel.

Example:

% if ndims == 2:
    fout[0][1] += t_xx;     fout[1][1] += t_xy;
    fout[0][2] += t_xy;     fout[1][2] += t_yy;
    fout[0][3] += u*t_xx + v*t_xy + ${-c['mu']*c['gamma']/c['Pr']}*T_x;
    fout[1][3] += u*t_xy + v*t_yy + ${-c['mu']*c['gamma']/c['Pr']}*T_y;
% endif

Loops

Mako loops can be used to facilitate PyFR-Mako kernel specification. Loops are opened with % for condition: and closed with % endfor. Note that such loops are unrolled when the PyFR-Mako kernel specification is interpreted at runtime, they are not embedded into the low-level kernel.

Example:

% for i in range(ndims):
    rhov[${i}] = s[${i + 1}];
    v[${i}] = invrho*rhov[${i}];
% endfor