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 physical-time
Controller available in PyFR 1.12.2:
StdNoneController Click to show
- class pyfr.integrators.std.controllers.StdNoneController(*args, **kwargs)[source]
- _accept_step(dt, idxcurr, err=None)
- _add(*args, subdims=None)
- _check_abort()
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_plugins()
- _get_reduction_kerns(*rs, **kwargs)
- _reject_step(dt, idxold, err=None)
- call_plugin_dt(dt)
- property cfgmeta
- collect_stats(stats)
- controller_name = 'none'
- property controller_needs_errest
- formulation = 'std'
- property grad_soln
- property nsteps
- run()
- property soln
- step(t, dt)
StdPIController Click to show
- class pyfr.integrators.std.controllers.StdPIController(*args, **kwargs)[source]
- _accept_step(dt, idxcurr, err=None)
- _add(*args, subdims=None)
- _check_abort()
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_plugins()
- _get_reduction_kerns(*rs, **kwargs)
- _reject_step(dt, idxold, err=None)
- call_plugin_dt(dt)
- property cfgmeta
- collect_stats(stats)
- controller_name = 'pi'
- property controller_needs_errest
- formulation = 'std'
- property grad_soln
- property nsteps
- run()
- property soln
- step(t, dt)
DualNoneController Click to show
- class pyfr.integrators.dual.phys.controllers.DualNoneController(*args, **kwargs)[source]
- _accept_step(idxcurr)
- _check_abort()
- _get_plugins()
- call_plugin_dt(dt)
- property cfgmeta
- collect_stats(stats)
- controller_name = 'none'
- formulation = 'dual'
- property grad_soln
- property nsteps
- property pseudostepinfo
- run()
- property soln
- step(t, dt)
- property system
Types of physical-time Controller are related via the following inheritance diagram:
There are two types of pseudo-time Controller available in PyFR 1.12.2:
DualNonePseudoController Click to show
- class pyfr.integrators.dual.pseudo.pseudocontrollers.DualNonePseudoController(*args, **kwargs)[source]
- _accumulate_source()
- _add(*args, subdims=None)
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_reduction_kerns(*rs, **kwargs)
- property _pseudo_stepper_regidx
- _resid(rcurr, rold, dt_fac)
- property _source_regidx
- property _stepper_regidx
- _update_pseudostepinfo(niters, resid)
- aux_nregs = 0
- convmon(i, minniters, dt_fac=1)
- finalise_pseudo_advance(currsoln)
- formulation = 'dual'
- pseudo_controller_name = 'none'
- pseudo_controller_needs_lerrest = False
DualPIPseudoController Click to show
- class pyfr.integrators.dual.pseudo.pseudocontrollers.DualPIPseudoController(*args, **kwargs)[source]
- _accumulate_source()
- _add(*args, subdims=None)
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_reduction_kerns(*rs, **kwargs)
- property _pseudo_stepper_regidx
- _resid(rcurr, rold, dt_fac)
- property _source_regidx
- property _stepper_regidx
- _update_pseudostepinfo(niters, resid)
- aux_nregs = 0
- convmon(i, minniters, dt_fac=1)
- finalise_pseudo_advance(currsoln)
- formulation = 'dual'
- pseudo_controller_name = 'local-pi'
- pseudo_controller_needs_lerrest = True
Types of pseudo-time Controller are related via the following inheritance diagram:
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 eight types of
Stepper available in PyFR 1.12.2:
StdEulerStepper Click to show
- class pyfr.integrators.std.steppers.StdEulerStepper(backend, systemcls, rallocs, mesh, initsoln, cfg)[source]
- _add(*args, subdims=None)
- _check_abort()
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_plugins()
- _get_reduction_kerns(*rs, **kwargs)
- property _stepper_nfevals
- advance_to(t)
- call_plugin_dt(dt)
- property cfgmeta
- collect_stats(stats)
- property controller_needs_errest
- formulation = 'std'
- property grad_soln
- property nsteps
- run()
- property soln
- stepper_has_errest = False
- stepper_name = 'euler'
- stepper_nregs = 2
- stepper_order = 1
StdRK4Stepper Click to show
- class pyfr.integrators.std.steppers.StdRK4Stepper(backend, systemcls, rallocs, mesh, initsoln, cfg)[source]
- _add(*args, subdims=None)
- _check_abort()
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_plugins()
- _get_reduction_kerns(*rs, **kwargs)
- property _stepper_nfevals
- advance_to(t)
- call_plugin_dt(dt)
- property cfgmeta
- collect_stats(stats)
- property controller_needs_errest
- formulation = 'std'
- property grad_soln
- property nsteps
- run()
- property soln
- stepper_has_errest = False
- stepper_name = 'rk4'
- stepper_nregs = 3
- stepper_order = 4
StdRK34Stepper Click to show
- class pyfr.integrators.std.steppers.StdRK34Stepper(*args, **kwargs)[source]
- _add(*args, subdims=None)
- _check_abort()
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_plugins()
- _get_reduction_kerns(*rs, **kwargs)
- _get_rkvdh2_kerns(stage, r1, r2, rold=None, rerr=None)
- property _stepper_nfevals
- 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)
- property cfgmeta
- collect_stats(stats)
- property controller_needs_errest
- formulation = 'std'
- property grad_soln
- property nsteps
- run()
- property soln
- step(t, dt)
- property stepper_has_errest
- stepper_name = 'rk34'
- property stepper_nregs
- stepper_order = 3
StdRK45Stepper Click to show
- class pyfr.integrators.std.steppers.StdRK45Stepper(*args, **kwargs)[source]
- _add(*args, subdims=None)
- _check_abort()
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_plugins()
- _get_reduction_kerns(*rs, **kwargs)
- _get_rkvdh2_kerns(stage, r1, r2, rold=None, rerr=None)
- property _stepper_nfevals
- 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)
- property cfgmeta
- collect_stats(stats)
- property controller_needs_errest
- formulation = 'std'
- property grad_soln
- property nsteps
- run()
- property soln
- step(t, dt)
- property stepper_has_errest
- stepper_name = 'rk45'
- property stepper_nregs
- stepper_order = 4
StdTVDRK3Stepper Click to show
- class pyfr.integrators.std.steppers.StdTVDRK3Stepper(backend, systemcls, rallocs, mesh, initsoln, cfg)[source]
- _add(*args, subdims=None)
- _check_abort()
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_plugins()
- _get_reduction_kerns(*rs, **kwargs)
- property _stepper_nfevals
- advance_to(t)
- call_plugin_dt(dt)
- property cfgmeta
- collect_stats(stats)
- property controller_needs_errest
- formulation = 'std'
- property grad_soln
- property nsteps
- run()
- property soln
- stepper_has_errest = False
- stepper_name = 'tvd-rk3'
- stepper_nregs = 3
- stepper_order = 3
DualBDF2Stepper Click to show
- class pyfr.integrators.dual.phys.steppers.DualBDF2Stepper(backend, systemcls, rallocs, mesh, initsoln, cfg)[source]
- _check_abort()
- _get_plugins()
- advance_to(t)
- call_plugin_dt(dt)
- property cfgmeta
- collect_stats(stats)
- formulation = 'dual'
- property grad_soln
- property nsteps
- property pseudostepinfo
- run()
- property soln
- step(t, dt)
- stepper_coeffs = [-1.5, 2.0, -0.5]
- stepper_name = 'bdf2'
- stepper_order = 2
- property system
DualBDF3Stepper Click to show
- class pyfr.integrators.dual.phys.steppers.DualBDF3Stepper(backend, systemcls, rallocs, mesh, initsoln, cfg)[source]
- _check_abort()
- _get_plugins()
- advance_to(t)
- call_plugin_dt(dt)
- property cfgmeta
- collect_stats(stats)
- formulation = 'dual'
- property grad_soln
- property nsteps
- property pseudostepinfo
- run()
- property soln
- step(t, dt)
- stepper_coeffs = [-1.8333333333333333, 3.0, -1.5, 0.3333333333333333]
- stepper_name = 'bdf3'
- stepper_order = 3
- property system
DualBackwardEulerStepper Click to show
- class pyfr.integrators.dual.phys.steppers.DualBackwardEulerStepper(backend, systemcls, rallocs, mesh, initsoln, cfg)[source]
- _check_abort()
- _get_plugins()
- advance_to(t)
- call_plugin_dt(dt)
- property cfgmeta
- collect_stats(stats)
- formulation = 'dual'
- property grad_soln
- property nsteps
- property pseudostepinfo
- run()
- property soln
- step(t, dt)
- stepper_coeffs = [-1.0, 1.0]
- stepper_name = 'backward-euler'
- stepper_order = 1
- property system
Types of Stepper are related via the following inheritance diagram:
PseudoStepper
A PseudoStepper acts to advance the simulation by a single pseudo-time-step. They are used to converge implicit Stepper time-steps via a dual time-stepping formulation. There are six types of PseudoStepper available in PyFR 1.12.2:
DualDenseRKPseudoStepper Click to show
- class pyfr.integrators.dual.pseudo.pseudosteppers.DualDenseRKPseudoStepper(*args, **kwargs)[source]
- _accumulate_source()
- _add(*args, subdims=None)
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_reduction_kerns(*rs, **kwargs)
- property _pseudo_stepper_regidx
- _rhs_with_dts(t, uin, fout)
- property _source_regidx
- property _stepper_regidx
- aux_nregs = 0
- collect_stats(stats)
- finalise_pseudo_advance(currsoln)
- formulation = 'dual'
- property ntotiters
- property pseudo_stepper_nfevals
DualRK4PseudoStepper Click to show
- class pyfr.integrators.dual.pseudo.pseudosteppers.DualRK4PseudoStepper(backend, systemcls, rallocs, mesh, initsoln, cfg, stepper_coeffs, dt)[source]
- _accumulate_source()
- _add(*args, subdims=None)
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_reduction_kerns(*rs, **kwargs)
- property _pseudo_stepper_regidx
- _rhs_with_dts(t, uin, fout)
- property _source_regidx
- property _stepper_regidx
- aux_nregs = 0
- collect_stats(stats)
- finalise_pseudo_advance(currsoln)
- formulation = 'dual'
- property ntotiters
- pseudo_stepper_has_lerrest = False
- pseudo_stepper_name = 'rk4'
- property pseudo_stepper_nfevals
- pseudo_stepper_nregs = 3
- pseudo_stepper_order = 4
DualTVDRK3PseudoStepper Click to show
- class pyfr.integrators.dual.pseudo.pseudosteppers.DualTVDRK3PseudoStepper(backend, systemcls, rallocs, mesh, initsoln, cfg, stepper_coeffs, dt)[source]
- _accumulate_source()
- _add(*args, subdims=None)
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_reduction_kerns(*rs, **kwargs)
- property _pseudo_stepper_regidx
- _rhs_with_dts(t, uin, fout)
- property _source_regidx
- property _stepper_regidx
- aux_nregs = 0
- collect_stats(stats)
- finalise_pseudo_advance(currsoln)
- formulation = 'dual'
- property ntotiters
- pseudo_stepper_has_lerrest = False
- pseudo_stepper_name = 'tvd-rk3'
- property pseudo_stepper_nfevals
- pseudo_stepper_nregs = 3
- pseudo_stepper_order = 3
DualEulerPseudoStepper Click to show
- class pyfr.integrators.dual.pseudo.pseudosteppers.DualEulerPseudoStepper(backend, systemcls, rallocs, mesh, initsoln, cfg, stepper_coeffs, dt)[source]
- _accumulate_source()
- _add(*args, subdims=None)
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_reduction_kerns(*rs, **kwargs)
- property _pseudo_stepper_regidx
- _rhs_with_dts(t, uin, fout)
- property _source_regidx
- property _stepper_regidx
- aux_nregs = 0
- collect_stats(stats)
- finalise_pseudo_advance(currsoln)
- formulation = 'dual'
- property ntotiters
- pseudo_stepper_has_lerrest = False
- pseudo_stepper_name = 'euler'
- property pseudo_stepper_nfevals
- pseudo_stepper_nregs = 2
- pseudo_stepper_order = 1
DualRK34PseudoStepper Click to show
- class pyfr.integrators.dual.pseudo.pseudosteppers.DualRK34PseudoStepper(*args, **kwargs)[source]
- _accumulate_source()
- _add(*args, subdims=None)
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_reduction_kerns(*rs, **kwargs)
- _get_rkvdh2pseudo_kerns(stage, r1, r2, rold, rerr=None)
- property _pseudo_stepper_regidx
- _rhs_with_dts(t, uin, fout)
- property _source_regidx
- property _stepper_regidx
- a = [0.32416573882874605, 0.5570978645055429, -0.08605491431272755]
- aux_nregs = 0
- b = [0.10407986927510238, 0.6019391368822611, 2.9750900268840206, -2.681109033041384]
- bhat = [0.3406814840808433, 0.09091523008632837, 2.866496742725443, -2.298093456892615]
- collect_stats(stats)
- finalise_pseudo_advance(currsoln)
- formulation = 'dual'
- property ntotiters
- property pseudo_stepper_has_lerrest
- pseudo_stepper_name = 'rk34'
- property pseudo_stepper_nfevals
- property pseudo_stepper_nregs
- pseudo_stepper_order = 3
- step(t)
DualRK45PseudoStepper Click to show
- class pyfr.integrators.dual.pseudo.pseudosteppers.DualRK45PseudoStepper(*args, **kwargs)[source]
- _accumulate_source()
- _add(*args, subdims=None)
- _get_axnpby_kerns(*rs, subdims=None)
- _get_gndofs()
- _get_reduction_kerns(*rs, **kwargs)
- _get_rkvdh2pseudo_kerns(stage, r1, r2, rold, rerr=None)
- property _pseudo_stepper_regidx
- _rhs_with_dts(t, uin, fout)
- property _source_regidx
- property _stepper_regidx
- a = [0.22502245872571303, 0.5440433129514047, 0.14456824349399464, 0.7866643421983568]
- aux_nregs = 0
- b = [0.05122930664033915, 0.3809548257264019, -0.3733525963923833, 0.5925012850263623, 0.34866717899927996]
- bhat = [0.13721732210321927, 0.19188076232938728, -0.2292067211595315, 0.6242946765438954, 0.27581396018302956]
- collect_stats(stats)
- finalise_pseudo_advance(currsoln)
- formulation = 'dual'
- property ntotiters
- property pseudo_stepper_has_lerrest
- pseudo_stepper_name = 'rk45'
- property pseudo_stepper_nfevals
- property pseudo_stepper_nregs
- pseudo_stepper_order = 4
- step(t)
Note that DualDenseRKPseudoStepper includes families of PseudoStepper whose coefficients are read from .txt files named thus:
{scheme name}-s{stage count}-p{temporal order}-sp{optimal spatial polynomial order}.txt
Types of PseudoStepper are related via the following inheritance diagram:
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.12.2 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 four types of System available
in PyFR 1.12.2:
ACEulerSystem Click to show
- class pyfr.solvers.aceuler.system.ACEulerSystem(backend, rallocs, mesh, initsoln, nregs, cfg)[source]
- _gen_kernels(eles, iint, mpiint, bcint)
- _gen_queues()
- _load_bc_inters(rallocs, mesh, elemap)
- _load_eles(rallocs, mesh, initsoln, nregs, nonce)
- _load_int_inters(rallocs, mesh, elemap)
- _load_mpi_inters(rallocs, mesh, elemap)
- _nonce_seq = count(0)
- _nqueues = 2
- bbcinterscls
alias of
pyfr.solvers.aceuler.inters.ACEulerBaseBCInters
- compute_grads(t, uinbank)
- ele_scal_upts(idx)
- elementscls
- filt(uinoutbank)
- intinterscls
- mpiinterscls
- name = 'ac-euler'
- rhs(t, uinbank, foutbank)
ACNavierStokesSystem Click to show
alias of
pyfr.solvers.acnavstokes.inters.ACNavierStokesBaseBCInters
EulerSystem Click to show
- class pyfr.solvers.euler.system.EulerSystem(backend, rallocs, mesh, initsoln, nregs, cfg)[source]
- _gen_kernels(eles, iint, mpiint, bcint)
- _gen_queues()
- _load_bc_inters(rallocs, mesh, elemap)
- _load_eles(rallocs, mesh, initsoln, nregs, 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
- compute_grads(t, uinbank)
- ele_scal_upts(idx)
- elementscls
- filt(uinoutbank)
- intinterscls
- mpiinterscls
- name = 'euler'
- rhs(t, uinbank, foutbank)
NavierStokesSystem Click to show
alias of
pyfr.solvers.navstokes.inters.NavierStokesBaseBCInters
Types of System are related via the following inheritance diagram:
Elements
An Elements holds information/data for a group of elements. There are four types of Elements available in PyFR 1.12.2:
ACEulerElements Click to show
- class pyfr.solvers.aceuler.elements.ACEulerElements(basiscls, eles, cfg)[source]
- _gen_pnorm_fpts()
- _mag_pnorm_fpts = None
- property _mesh_regions
- _norm_pnorm_fpts = None
- _ploc_in_src_exprs = None
- property _scratch_bufs
- _slice_mat(mat, region, ra=None, rb=None)
- _smats_djacs_mpts = None
- _soln_in_src_exprs = None
- _src_exprs = None
- _srtd_face_fpts = None
- static con_to_pri(convs, cfg)
- convarmap = {2: ['p', 'u', 'v'], 3: ['p', 'u', 'v', 'w']}
- dualcoeffs = {2: ['u', 'v'], 3: ['u', 'v', 'w']}
- formulations = ['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, side=None)
- ploc_at_np(name)
- plocfpts = None
- static pri_to_con(pris, cfg)
- privarmap = {2: ['p', 'u', 'v'], 3: ['p', 'u', 'v', 'w']}
- qpts = None
- rcpdjac_at(name, side=None)
- rcpdjac_at_np(name)
- set_ics_from_cfg()
- set_ics_from_soln(solnmat, solncfg)
- sliceat()
- smat_at(name, side=None)
- smat_at_np(name)
- upts = None
- visvarmap = {2: [('velocity', ['u', 'v']), ('pressure', ['p'])], 3: [('velocity', ['u', 'v', 'w']), ('pressure', ['p'])]}
ACNavierStokesElements Click to show
EulerElements Click to show
- class pyfr.solvers.euler.elements.EulerElements(basiscls, eles, cfg)[source]
- _gen_pnorm_fpts()
- _mag_pnorm_fpts = None
- property _mesh_regions
- _norm_pnorm_fpts = None
- _ploc_in_src_exprs = None
- property _scratch_bufs
- _slice_mat(mat, region, ra=None, rb=None)
- _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, side=None)
- ploc_at_np(name)
- plocfpts = None
- static pri_to_con(pris, cfg)
- privarmap = {2: ['rho', 'u', 'v', 'p'], 3: ['rho', 'u', 'v', 'w', 'p']}
- qpts = None
- rcpdjac_at(name, side=None)
- rcpdjac_at_np(name)
- set_ics_from_cfg()
- set_ics_from_soln(solnmat, solncfg)
- sliceat()
- smat_at(name, side=None)
- smat_at_np(name)
- upts = None
- visvarmap = {2: [('density', ['rho']), ('velocity', ['u', 'v']), ('pressure', ['p'])], 3: [('density', ['rho']), ('velocity', ['u', 'v', 'w']), ('pressure', ['p'])]}
NavierStokesElements Click to show
Types of Elements are related via the following inheritance diagram:
Interfaces
An Interfaces holds information/data for a group of interfaces. There are eight types of (non-boundary) Interfaces available in PyFR 1.12.2:
ACEulerIntInters Click to show
- class pyfr.solvers.aceuler.inters.ACEulerIntInters(*args, **kwargs)[source]
- _const_mat(inter, meth)
- _gen_perm(lhs, rhs)
- _get_perm_for_view(inter, meth)
- _scal_view(inter, meth)
- _scal_xchg_view(inter, meth)
- _set_external(name, spec, value=None)
- _vect_view(inter, meth)
- _vect_xchg_view(inter, meth)
- _view(inter, meth, vshape=())
- _xchg_view(inter, meth, vshape=())
- prepare(t)
ACEulerMPIInters Click to show
- class pyfr.solvers.aceuler.inters.ACEulerMPIInters(*args, **kwargs)[source]
- MPI_TAG = 2314
- _const_mat(inter, meth)
- _get_perm_for_view(inter, meth)
- _scal_view(inter, meth)
- _scal_xchg_view(inter, meth)
- _set_external(name, spec, value=None)
- _vect_view(inter, meth)
- _vect_xchg_view(inter, meth)
- _view(inter, meth, vshape=())
- _xchg_view(inter, meth, vshape=())
- prepare(t)
ACNavierStokesIntInters Click to show
ACNavierStokesMPIInters Click to show
EulerIntInters Click to show
- class pyfr.solvers.euler.inters.EulerIntInters(*args, **kwargs)[source]
- _const_mat(inter, meth)
- _gen_perm(lhs, rhs)
- _get_perm_for_view(inter, meth)
- _scal_view(inter, meth)
- _scal_xchg_view(inter, meth)
- _set_external(name, spec, value=None)
- _vect_view(inter, meth)
- _vect_xchg_view(inter, meth)
- _view(inter, meth, vshape=())
- _xchg_view(inter, meth, vshape=())
- prepare(t)
EulerMPIInters Click to show
- class pyfr.solvers.euler.inters.EulerMPIInters(*args, **kwargs)[source]
- MPI_TAG = 2314
- _const_mat(inter, meth)
- _get_perm_for_view(inter, meth)
- _scal_view(inter, meth)
- _scal_xchg_view(inter, meth)
- _set_external(name, spec, value=None)
- _vect_view(inter, meth)
- _vect_xchg_view(inter, meth)
- _view(inter, meth, vshape=())
- _xchg_view(inter, meth, vshape=())
- prepare(t)
NavierStokesIntInters Click to show
NavierStokesMPIInters Click to show
Types of (non-boundary) Interfaces are related via the following inheritance diagram:
Backend
A Backend holds information/data for a backend. There are four types of Backend available in PyFR 1.12.2:
CUDABackend Click to show
- class pyfr.backends.cuda.base.CUDABackend(cfg)[source]
-
- alias(obj, aobj)
- blocks = False
- 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_slice(mat, ra, rb, ca, cb)
- 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={})
HIPBackend Click to show
- class pyfr.backends.hip.base.HIPBackend(cfg)[source]
-
- alias(obj, aobj)
- blocks = False
- 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_slice(mat, ra, rb, ca, cb)
- name = 'hip'
- 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={})
OpenCLBackend Click to show
- class pyfr.backends.opencl.base.OpenCLBackend(cfg)[source]
-
- alias(obj, aobj)
- blocks = False
- 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_slice(mat, ra, rb, ca, cb)
- 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={})
OpenMPBackend Click to show
- class pyfr.backends.openmp.base.OpenMPBackend(cfg)[source]
-
- alias(obj, aobj)
- blocks = True
- 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_slice(mat, ra, rb, ca, cb)
- 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:
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.12.2:
CUDAPointwiseKernelProvider Click to show
HIPPointwiseKernelProvider Click to show
OpenCLPointwiseKernelProvider Click to show
- class pyfr.backends.opencl.provider.OpenCLPointwiseKernelProvider(backend)[source]
- _build_arglst(dims, argn, argt, argdict)
- _build_kernel(name, src, argtypes)
- _render_kernel(name, mod, extrns, tplargs)
- kernel_generator_cls
alias of
pyfr.backends.opencl.generator.OpenCLKernelGenerator
- register(mod)
OpenMPPointwiseKernelProvider Click to show
- class pyfr.backends.openmp.provider.OpenMPPointwiseKernelProvider(backend)[source]
- _build_arglst(dims, argn, argt, argdict)
- _build_kernel(name, src, argtypes, restype=None)
- _render_kernel(name, mod, extrns, 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:
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.12.2:
CUDAKernelGenerator Click to show
HIPKernelGenerator Click to show
OpenCLKernelGenerator Click to show
OpenMPKernelGenerator Click to show
Types of Kernel Generator are related via the following inheritance diagram:
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
kernel-name
— name of kernelstring
data-dimensionality
— dimensionality of dataint
argument-name
— name of argumentstring
argument-intent
— intent of argumentin
|out
|inout
argument-attribute
— attribute of argumentmpi
|scalar
|view
argument-data-type
— data type of argumentstring
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
macro-name
— name of macrostring
parameter-name
— name of parameterstring
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
macro-name
— name of the macrostring
parameter-name
— name of parameterstring
Syntax
Basic Functionality
Basic functionality can be expressed using a restricted subset of the C programming language. Specifically, use of the following is allowed:
+,-,*,/
— basic arithmeticsin, cos, tan
— basic trigonometric functionsexp
— exponentialpow
— powerfabs
— absolute valueoutput = ( condition ? satisfied : unsatisfied )
— ternary ifmin
— minimummax
— 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