Provided by: pyfr_1.5.0-3_all bug

NAME

       pyfr - PyFR Documentation

       Contents:

HOME

   Overview
   What is PyFR?
       PyFR  is an open-source Python based framework for solving advection-diffusion type problems on streaming
       architectures using the Flux Reconstruction approach of Huynh. The framework is designed to solve a range
       of governing systems on mixed unstructured grids containing various element types. It is also designed to
       target a range of hardware platforms via use of an in-built domain specific  language  derived  from  the
       Mako templating engine. The current release (PyFR 1.5.0) has the following capabilities:

       • Governing Equations - Euler, Navier Stokes

       • Dimensionality - 2D, 3D

       • Element Types - Triangles, Quadrilaterals, Hexahedra, Prisms, Tetrahedra, Pyramids

       • Platforms - CPU Clusters, Nvidia GPU Clusters, AMD GPU Clusters, Intel Xeon Phi Clusters

       • Spatial Discretisation - High-Order Flux Reconstruction

       • Temporal Discretisation - Explicit and Implicit (via Dual Time-Stepping)

       • Precision - Single, Double

       • Mesh Files Imported - Gmsh (.msh), CGNS (.cgns)

       • Solution Files Exported - Unstructured VTK (.vtu, .pvtu)

   How do I Cite PyFR?
       To cite PyFR, please reference the following paper:

       • PyFR: An Open Source Framework for Solving Advection-Diffusion Type Problems on Streaming Architectures
         using the Flux Reconstruction Approach. F. D. Witherden, A. M.  Farrington,  P.  E.  Vincent.  Computer
         Physics Communications, Volume 185, Pages 3028-3040, 2014.

   Who is Funding PyFR?
       Development  of PyFR is supported by the Engineering and Physical Sciences Research Council, Innovate UK,
       the European Commission, BAE Systems, Airbus, and the Air Force Office of Scientific  Research.   We  are
       also grateful for hardware donations from Nvidia, Intel, and AMD.

THEORY

   Flux Reconstruction
   Overview
       High-order  numerical methods for unstructured grids combine the superior accuracy of high-order spectral
       or finite difference methods with the geometrical  flexibility  of  low-order  finite  volume  or  finite
       element   schemes.  The  Flux  Reconstruction  (FR)  approach  unifies  various  high-order  schemes  for
       unstructured grids within a single framework.  Additionally,  the  FR  approach  exhibits  a  significant
       degree  of  element locality, and is thus able to run efficiently on modern streaming architectures, such
       as Graphical Processing Units (GPUs). The aforementioned properties of FR  mean  it  offers  a  promising
       route  to performing affordable, and hence industrially relevant, scale-resolving simulations of hitherto
       intractable unsteady flows (involving separation, acoustics  etc.)  within  the  vicinity  of  real-world
       engineering geometries. An detailed overview of the FR approach is given in:

       • A  Flux  Reconstruction  Approach to High-Order Schemes Including Discontinuous Galerkin Methods. H. T.
         Huynh. AIAA Paper 2007-4079

   Linear Stability
       The linear stability of an FR schemes depends on the form of the correction  function.  Linear  stability
       issues are discussed in:

       • A  New Class of High-Order Energy Stable Flux Reconstruction Schemes.  P. E. Vincent, P. Castonguay, A.
         Jameson. Journal of Scientific Computing, Volume 47, Number 1, Pages 50-72, 2011Insights from von Neumann Analysis of  High-Order  Flux  Reconstruction  Schemes.  P.  E.  Vincent,  P.
         Castonguay, A. Jameson. Journal of Computational Physics, Volume 230, Issue 22, Pages 8134-8154, 2011A  New  Class  of  High-Order  Energy  Stable  Flux  Reconstruction Schemes for Triangular Elements. P.
         Castonguay, P. E. Vincent, A. Jameson. Journal of Scientific Computing,  Volume  51,  Number  1,  Pages
         224-256, 2012Energy  Stable  Flux  Reconstruction  Schemes  for  Advection-Diffusion Problems.  P. Castonguay, D. M.
         Williams, P. E. Vincent, A. Jameson. Computer Methods in Applied Mechanics and Engineering, Volume 267,
         Pages 400-417, 2013Energy  Stable  Flux  Reconstruction  Schemes  for  Advection-Diffusion  Problems  on  Triangles. D. M.
         Williams, P. Castonguay, P. E. Vincent, A. Jameson.  Journal  of  Computational  Physics,  Volume  250,
         Pages 53-76, 2013Energy  Stable  Flux  Reconstruction  Schemes  for  Advection-Diffusion  Problems  on Tetrahedra. D. M.
         Williams, A. Jameson. Journal of Scientific Computing, Volume 59, Pages 721-759, 2014An Extended Range of Stable-Symmetric-Conservative Flux  Reconstruction  Correction  Functions.  P.  E.
         Vincent,  A.  M.  Farrington,  F.  D.  Witherden, A. Jameson. Computer Methods in Applied Mechanics and
         Engineering, Volume 296, Pages 248-272, 2015

   Non-Linear Stability
       The non-linear stability of an FR schemes depends on the location  of  the  solution  points.  Non-linear
       stability issues are discussed in:

       • On  the  Non-Linear Stability of Flux Reconstruction Schemes. A. Jameson, P. E. Vincent, P. Castonguay.
         Journal of Scientific Computing, Volume 50, Number 2, Pages 434-445, 2012An Analysis of Solution Point Coordinates for Flux Reconstruction Schemes on Triangular Elements. F. D.
         Witherden, P. E. Vincent. Journal of Scientific Computing, Volume 61, Pages 398-423, 2014

USER GUIDE

   Getting Started
   Downloading the Source
       PyFR can be obtained here.

   Dependencies
   Overview
       PyFR 1.5.0 has a hard dependency on Python 3.3+ and the following Python packages:

       1. gimmik >= 2.0

       2. h5py >= 2.6

       3. mako >= 1.0.0

       4. mpi4py >= 2.0

       5. numpy >= 1.8

       6. pytools >= 2016.2.1

       Note that due to a bug in numpy PyFR is not compatible with 32-bit Python distributions.

   CUDA Backend
       The CUDA backend targets NVIDIA GPUs with a compute capability of 2.0 or greater. The backend requires:

       1. CUDA >= 4.2

       2. pycuda >= 2015.1

   MIC Backend
       The MIC backend targets Intel Xeon Phi co-processors. The backend requires:

       1. ICC >= 14.0

       2. Intel MKL >= 11.1

       3. Intel MPSS >= 3.3

       4. pymic >= 0.7 (post commit 4d8a2da)

   OpenCL Backend
       The  OpenCL  backend  targets  a  range  of  accelerators including GPUs from AMD and NVIDIA. The backend
       requires:

       1. OpenCL

       2. pyopencl >= 2015.2.4

       3. clBLAS

   OpenMP Backend
       The OpenMP backend targets multi-core CPUs. The backend requires:

       1. GCC >= 4.9

       2. A BLAS library compiled as a shared library (e.g. OpenBLAS)

   Running in Parallel
       To partition meshes for running  in  parallel  it  is  also  necessary  to  have  one  of  the  following
       partitioners installed:

       1. metis >= 5.0

       2. scotch >= 6.0

   Importing CGNS Meshes
       To import CGNS meshes it is necessary to have the following installed:

       1. CGNS >= 3.3 (develop branch post commit e0faea6)

   Installation
       Before  running  PyFR  1.5.0  it  is  first  necessary  to either install the software using the provided
       setup.py installer or add the root PyFR directory to PYTHONPATH using:

          user@computer ~/PyFR$ export PYTHONPATH=.:$PYTHONPATH

       To manage installation of Python dependencies we strongly recommend using pip and virtualenv.

   Running PyFR
   Overview
       PyFR 1.5.0 uses three distinct file formats:

       1. .ini --- configuration file

       2. .pyfrm --- mesh file

       3. .pyfrs --- solution file

       The following commands are available from the pyfr program:

       1. pyfr import --- convert a Gmsh .msh file or CGNS .cgns file into a PyFR .pyfrm file.

          Example:

             pyfr import mesh.msh mesh.pyfrm

       2. pyfr partition --- partition an existing mesh and associated solution files.

          Example:

             pyfr partition 2 mesh.pyfrm solution.pyfrs .

       3. pyfr run --- start a new PyFR simulation. Example:

             pyfr run mesh.pyfrm configuration.ini

       4. pyfr restart --- restart a PyFR simulation from an existing solution file. Example:

             pyfr restart mesh.pyfrm solution.pyfrs

       5. pyfr export --- convert a PyFR .pyfrs file into an unstructured VTK .vtu or .pvtu file. Example:

             pyfr export mesh.pyfrm solution.pyfrs solution.vtu

   Running in Parallel
       pyfr can be run in parallel. To do so prefix pyfr with mpirun -n <cores/devices>. Note that the mesh must
       be pre-partitioned, and the number of cores or devices must be equal to the number of partitions.

   Configuration File (.ini)
   Overview
       The  .ini  configuration  file parameterises the simulation. It is written in the INI format.  Parameters
       are grouped into sections. The roles of each section and their associated parameters are described below.

   [backend]
       Parameterises the backend with

       1. precision --- number precision:
             single | double

       2. rank-allocator --- MPI rank allocator:
             linear | random

       Example:

          [backend]
          precision = double
          rank-allocator = linear

   [backend-cuda]
       Parameterises the CUDA backend with

       1. device-id --- method for selecting which device(s) to run on:
             int | round-robin | local-rank

       2. gimmik-max-nnz --- cutoff for GiMMiK in terms of the number of non-zero entires in a constant matrix:
             int

       3. mpi-type --- type of MPI library that is being used:
             standard | cuda-aware

       4. block-1d --- block size for one dimensional pointwise kernels:
             int

       5. block-2d --- block size for two dimensional pointwise kernels:
             int, int

       Example:

          [backend-cuda]
          device-id = round-robin
          gimmik-max-nnz = 512
          mpi-type = standard
          block-1d = 64
          block-2d = 128, 2

   [backend-mic]
       Parameterises the MIC backend with

       1. device-id --- for selecting which device(s) to run on:
             int | local-rank

       2. mkl-root --- path to MKL root directory:
             string

   [backend-opencl]
       Parameterises the OpenCL backend with

       1. platform-id --- for selecting platform id:
             int | string

       2. device-type --- for selecting what type of device(s) to run on:
             all | cpu | gpu | accelerator

       3. device-id --- for selecting which device(s) to run on:
             int | string | local-rank

       4. gimmik-max-nnz --- cutoff for GiMMiK in terms of the number of non-zero entires in a constant matrix:
             int

       5. local-size-1d --- local work size for one dimensional pointwise kernels:
             int

       6. local-size-2d --- local work size for two dimensional pointwise kernels:
             int, int

       Example:

          [backend-opencl]
          platform-id = 0
          device-type = gpu
          device-id = local-rank
          gimmik-max-nnz = 512
          local-size-1d = 16
          local-size-2d = 128, 1

   [backend-openmp]
       Parameterises the OpenMP backend with

       1. cc --- C compiler:
             string

       2. cflags --- additional C compiler flags:
             string

       3. cblas --- path to shared C BLAS library:
             string

       4. cblas-type --- type of BLAS library:
             serial | parallel

       Example:

          [backend-openmp]
          cc = gcc
          cblas= example/path/libBLAS.dylib
          cblas-type = parallel

   [constants]
       Sets constants used in the simulation with

       1. gamma --- ratio of specific heats:
             float

       2. mu --- dynamic viscosity:
             float

       3. Pr --- Prandtl number:
             float

       4. cpTref --- product of specific heat at constant pressure and reference  temperature  for  Sutherland's
          Law:

          float

       5. cpTs  ---  product  of  specific heat at constant pressure and Sutherland temperature for Sutherland's
          Law:

          float

       Example:

          [constants]
          gamma = 1.4
          mu = 0.001
          Pr = 0.72

   [solver]
       Parameterises the solver with

       1. system --- governing system:
             euler | navier-stokes

       2. order --- order of polynomial solution basis:
             int

       3. anti-alias --- type of anti-aliasing:
             flux | surf-flux | div-flux | flux, surf-flux |  flux,  div-flux  |  surf-flux,  div-flux  |  flux,
             surf-flux, div-flux

       4. viscosity-correction --- viscosity correction:
             none | sutherland

       5. shock-capturing --- shock capturing scheme:
             none | artificial-viscosity

       Example:

          [solver]
          system = navier-stokes
          order = 3
          anti-alias = flux
          viscosity-correction = none
          shock-capturing = artificial-viscosity

   [solver-time-integrator]
       Parameterises the time-integration scheme used by the solver with

       1. formulation --- formulation:
             std | dual

             where

             std requires

                 • scheme --- time-integration scheme
                       euler | rk34 | rk4 | rk45 | tvd-rk3tstart --- initial time
                       floattend --- final time
                       floatdt --- time-step
                       floatcontroller --- time-step controller
                       none | pi

                       where

                       pi only works with rk34 and rk45 and requires

                           • atol --- absolute error tolerance
                                floatrtol --- relative error tolerance
                                floaterrest-norm --- norm to use for estimating the error
                                uniform | l2safety-fact --- safety factor for step size adjustment (suitable range 0.80-0.95)
                                floatmin-fact  ---  minimum  factor  that  the  time-step  can change between iterations
                             (suitable range 0.1-0.5)
                                floatmax-fact --- maximum factor  that  the  time-step  can  change  between  iterations
                             (suitable range 2.0-6.0)
                                float

             dual requires

                 • scheme --- time-integration scheme
                       backward-euler | bdf2 | bdf3pseudo-scheme --- pseudo-time-integration scheme
                       euler | tvd-rk3 | rk4tstart --- initial time
                       floattend --- final time
                       floatdt --- time-step
                       floatpseudo-dt --- pseudo-time-step
                       floatcontroller --- pseudo-time-step controller
                       none

                       where

                       none requires

                           • pseudo-niters-max --- minimum number of iterations
                                intpseudo-niters-min --- maximum number of iterations
                                intpseudo-aresid --- absolute residual tolerance
                                floatpseudo-rresid --- relative residual tolerance
                                float

       Example:

          [solver-time-integrator]
          formulation = std
          scheme = rk45
          controller = pi
          tstart = 0.0
          tend = 10.0
          dt = 0.001
          atol = 0.00001
          rtol = 0.00001
          errest-norm = l2
          safety-fact = 0.9
          min-fact = 0.3
          max-fact = 2.5

   [solver-interfaces]
       Parameterises the interfaces with

       1. riemann-solver --- type of Riemann solver:
             rusanov | hll | hllc | roe | roem

       2. ldg-beta --- beta parameter used for LDG:
             float

       3. ldg-tau --- tau parameter used for LDG:
             float

       Example:

          [solver-interfaces]
          riemann-solver = rusanov
          ldg-beta = 0.5
          ldg-tau = 0.1

   [solver-interfaces-line]
       Parameterises the line interfaces with

       1. flux-pts --- location of the flux points on a line interface:
             gauss-legendre | gauss-legendre-lobatto

       2. quad-deg --- degree of quadrature rule for anti-aliasing on a line interface:
             int

       3. quad-pts --- name of quadrature rule for anti-aliasing on a line interface:
             gauss-legendre | gauss-legendre-lobatto

       Example:

          [solver-interfaces-line]
          flux-pts = gauss-legendre
          quad-deg = 10
          quad-pts = gauss-legendre

   [solver-interfaces-tri]
       Parameterises the triangular interfaces with

       1. flux-pts --- location of the flux points on a triangular interface:
             williams-shunn

       2. quad-deg --- degree of quadrature rule for anti-aliasing on a triangular interface:
             int

       3. quad-pts --- name of quadrature rule for anti-aliasing on a triangular interface:
             williams-shunn | witherden-vincent

       Example:

          [solver-interfaces-tri]
          flux-pts = williams-shunn
          quad-deg = 10
          quad-pts = williams-shunn

   [solver-interfaces-quad]
       Parameterises the quadrilateral interfaces with

       1. flux-pts --- location of the flux points on a quadrilateral interface:
             gauss-legendre | gauss-legendre-lobatto

       2. quad-deg --- degree of quadrature rule for anti-aliasing on a quadrilateral interface:
             int

       3. quad-pts --- name of quadrature rule for anti-aliasing on a quadrilateral interface:
             gauss-legendre | gauss-legendre-lobatto | witherden-vincent

       Example:

          [solver-interfaces-quad]
          flux-pts = gauss-legendre
          quad-deg = 10
          quad-pts = gauss-legendre

   [solver-elements-tri]
       Parameterises the triangular elements with

       1. soln-pts --- location of the solution points in a triangular element:
             williams-shunn

       2. quad-deg --- degree of quadrature rule for anti-aliasing in a triangular element:
             int

       3. quad-pts --- name of quadrature rule for anti-aliasing in a triangular element:
             williams-shunn | witherden-vincent

       Example:

          [solver-elements-tri]
          soln-pts = williams-shunn
          quad-deg = 10
          quad-pts = williams-shunn

   [solver-elements-quad]
       Parameterises the quadrilateral elements with

       1. soln-pts --- location of the solution points in a quadrilateral element:
             gauss-legendre | gauss-legendre-lobatto

       2. quad-deg --- degree of quadrature rule for anti-aliasing in a quadrilateral element:
             int

       3. quad-pts --- name of quadrature rule for anti-aliasing in a quadrilateral element:
             gauss-legendre | gauss-legendre-lobatto | witherden-vincent

       Example:

          [solver-elements-quad]
          soln-pts = gauss-legendre
          quad-deg = 10
          quad-pts = gauss-legendre

   [solver-elements-hex]
       Parameterises the hexahedral elements with

       1. soln-pts --- location of the solution points in a hexahedral element:
             gauss-legendre | gauss-legendre-lobatto

       2. quad-deg --- degree of quadrature rule for anti-aliasing in a hexahedral element:
             int

       3. quad-pts --- name of quadrature rule for anti-aliasing in a hexahedral element:
             gauss-legendre | gauss-legendre-lobatto | witherden-vincent

       Example:

          [solver-elements-hex]
          soln-pts = gauss-legendre
          quad-deg = 10
          quad-pts = gauss-legendre

   [solver-elements-tet]
       Parameterises the tetrahedral elements with

       1. soln-pts --- location of the solution points in a tetrahedral element:
             shunn-ham

       2. quad-deg --- degree of quadrature rule for anti-aliasing in a tetrahedral element:
             int

       3. quad-pts --- name of quadrature rule for anti-aliasing in a tetrahedral element:
             shunn-ham | witherden-vincent

       Example:

          [solver-elements-tet]
          soln-pts = shunn-ham
          quad-deg = 10
          quad-pts = shunn-ham

   [solver-elements-pri]
       Parameterises the prismatic elements with

       1. soln-pts --- location of the solution points in a prismatic element:
             williams-shunn~gauss-legendre | williams-shunn~gauss-legendre-lobatto

       2. quad-deg --- degree of quadrature rule for anti-aliasing in a prismatic element:
             int

       3. quad-pts --- name of quadrature rule for anti-aliasing in a prismatic element:
             williams-shunn~gauss-legendre | williams-shunn~gauss-legendre-lobatto | witherden-vincent

       Example:

          [solver-elements-pri]
          soln-pts = williams-shunn~gauss-legendre
          quad-deg = 10
          quad-pts = williams-shunn~gauss-legendre

   [solver-elements-pyr]
       Parameterises the pyramidal elements with

       1. soln-pts --- location of the solution points in a pyramidal element:
             gauss-legendre | gauss-legendre-lobatto

       2. quad-deg --- degree of quadrature rule for anti-aliasing in a pyramidal element:
             int

       3. quad-pts --- name of quadrature rule for anti-aliasing in a pyramidal element:
             witherden-vincent

       Example:

          [solver-elements-pyr]
          soln-pts = gauss-legendre
          quad-deg = 10
          quad-pts = witherden-vincent

   [solver-source-terms]
       Parameterises solution, space (x, y, [z]), and time (t) dependent source terms with

       1. rho --- density source term:
             string

       2. rhou --- x-momentum source term:
             string

       3. rhov --- y-momentum source term:
             string

       4. rhow --- z-momentum source term:
             string

       5. E --- energy source term:
             string

       Example:

          [solver-source-terms]
          rho = t
          rhou = x*y*sin(y)
          rhov = z*rho
          rhow = 1.0
          E = 1.0/(1.0+x)

   [solver-artificial-viscosity]
       Parameterises artificial viscosity for shock capturing with

       1. max-artvisc --- maximum artificial viscosity:
             float

       2. s0 --- sensor cut-off:
             float

       3. kappa --- sensor range:
             float

       Example:

          [solver-artificial-viscosity]
          max-artvisc = 0.01
          s0 = 0.01
          kappa = 5.0

   [soln-filter]
       Parameterises an exponential solution filter with

       1. nsteps --- apply filter every nsteps:
             int

       2. alpha --- strength of filter:
             float

       3. order --- order of filter:
             int

       4. cutoff --- cutoff frequency below which no filtering is applied:
             int

       Example:

          [soln-filter]
          nsteps = 10
          alpha = 36.0
          order = 16
          cutoff = 1

   [soln-plugin-writer]
       Periodically write the solution to disk in the pyfrs format.  Parameterised with

       1. dt-out --- write to disk every dt-out time units:
             float

       2. basedir --- relative path to directory where outputs will be written:
             string

       3. basename --- pattern of output names:
             string

       4. post-action --- command to execute after writing the file:
             string

       5. post-action-mode --- how the post-action command should be executed:
             blocking | non-blocking

       Example:

          [soln-plugin-writer]
          dt-out = 0.01
          basedir = .
          basename = files-{t:.2f}
          post-action = echo "Wrote file {soln} at time {t} for mesh {mesh}."
          post-action-mode = blocking

   [soln-plugin-fluidforce-name]
       Periodically  integrates the pressure and viscous stress on the boundary labelled name and writes out the
       resulting force vectors to a CSV file. Parameterised with

       1. nsteps --- integrate every nsteps:
             int

       2. file --- output file path; should the file already exist it will be appended to:
             string

       3. header --- if to output a header row or not:
             boolean

       Example:

          [soln-plugin-fluidforce-wing]
          nsteps = 10
          file = wing-forces.csv
          header = true

   [soln-plugin-nancheck]
       Periodically checks the solution for NaN values. Parameterised with

       1. nsteps --- check every nsteps:
             int

       Example:

          [soln-plugin-nancheck]
          nsteps = 10

   [soln-plugin-residual]
       Periodically calculates the residual and writes it out to a CSV file.  Parameterised with

       1. nsteps --- calculate every nsteps:
             int

       2. file --- output file path; should the file already exist it will be appended to:
             string

       3. header --- if to output a header row or not:
             boolean

       Example:

          [soln-plugin-residual]
          nsteps = 10
          file = residual.csv
          header = true

   [soln-plugin-dtstats]
       Write time-step statistics out to a CSV file. Parameterised with

       1. flushsteps --- flush to disk every flushsteps:
             int

       2. file --- output file path; should the file already exist it will be appended to:
             string

       3. header --- if to output a header row or not:
             boolean

       Example:

          [soln-plugin-dtstats]
          flushsteps = 100
          file = dtstats.csv
          header = true

   [soln-plugin-sampler]
       Periodically samples specific points in the volume and writes  them  out  to  a  CSV  file.   The  plugin
       actually  samples  the  solution  point  closest  to each sample point, hence a slight discrepancy in the
       output sampling locations is to be expected.  A nearest-neighbour search is used to  locate  the  closest
       solution   point   to   the  sample  point.   The  location  process  automatically  takes  advantage  of
       scipy.spatial.cKDTree where available.  Parameterised with

       1. nsteps --- sample every nsteps:
             int

       2. samp-pts --- list of points to sample:
             [(x, y), (x, y), ...] | [(x, y, z), (x, y, z), ...]

       3. format --- output variable format:
             primitive | conservative

       4. file --- output file path; should the file already exist it will be appended to:
             string

       5. header --- if to output a header row or not:
             boolean

       Example:

          [soln-plugin-sampler]
          nsteps = 10
          samp-pts = [(1.0, 0.7, 0.0), (1.0, 0.8, 0.0)]
          format = primative
          file = point-data.csv
          header = true

   [soln-plugin-tavg]
       Time average quantities. Parameterised with

       1. nsteps --- accumulate the average every nsteps time steps:
             int

       2. dt-out --- write to disk every dt-out time units:
             float

       3. basedir --- relative path to directory where outputs will be written:
             string

       4. basename --- pattern of output names:
             string

       5. avg-name --- expression as a function of the primitive variables, time (t), and space (x, y,  [z])  to
          time average; multiple expressions, each with their own name, may be specified:
             string

       Example:

          [soln-plugin-tavg]
          nsteps = 10
          dt-out = 2.0
          basedir = .
          basename = files-{t:06.2f}

          avg-p = p
          avg-p2 = p*p
          avg-vel = sqrt(u*u + v*v)

   [soln-bcs-name]
       Parameterises  constant,  or  if  available  space (x, y, [z]) and time (t) dependent, boundary condition
       labelled name in the .pyfrm file with

       1. type --- type of boundary condition:
             char-riem-inv | no-slp-adia-wall | no-slp-isot-wall | slp-adia-wall | sub-in-frv |  sub-in-ftpttang
             | sub-out-fp | sup-in-fa | sup-out-fn

             where

             char-riem-inv requires

                 • rho --- density
                       float | stringu --- x-velocity
                       float | stringv --- y-velocity
                       float | stringw --- z-velocity
                       float | stringp --- static pressure
                       float | string

             no-slp-isot-wall requires

                 • u --- x-velocity of wall
                       floatv --- y-velocity of wall
                       floatw --- z-velocity of wall
                       floatcpTw --- product of specific heat capacity at constant pressure and temperature of wall
                       float

             sub-in-frv requires

                 • rho --- density
                       float | stringu --- x-velocity
                       float | stringv --- y-velocity
                       float | stringw --- z-velocity
                       float | string

             sub-in-ftpttang requires

                 • pt --- total pressure
                       floatcpTt --- product of specific heat capacity at constant pressure and total temperature
                       floattheta  ---  azimuth  angle  (in  degrees) of inflow measured in the x-y plane relative to the
                   positive x-axis
                       floatphi --- inclination angle (in degrees) of inflow measured relative to the positive z-axis
                       float

             sub-out-fp requires

                 • p --- static pressure
                       float | string

             sup-in-fa requires

                 • rho --- density
                       float | stringu --- x-velocity
                       float | stringv --- y-velocity
                       float | stringw --- z-velocity
                       float | stringp --- static pressure
                       float | string

       Example:

          [soln-bcs-bcwallupper]
          type = no-slp-isot-wall
          cpTw = 10.0
          u = 1.0

   [soln-ics]
       Parameterises space (x, y, [z]) dependent initial conditions with

       1. rho --- initial density distribution:
             string

       2. u --- initial x-velocity distribution:
             string

       3. v --- initial y-velocity distribution:
             string

       4. w --- initial z-velocity distribution:
             string

       5. p --- initial static pressure distribution:
             string

       Example:

          [soln-ics]
          rho = 1.0
          u = x*y*sin(y)
          v = z
          w = 1.0
          p = 1.0/(1.0+x)

   Example --- 2D Couette Flow
       Proceed with the following steps to run a serial 2D Couette flow simulation on a mixed unstructured mesh:

       1. Create a working directory called couette_flow_2d/

       2. Copy the configuration file PyFR/examples/couette_flow_2d/couette_flow_2d.ini into couette_flow_2d/

       3. Copy the Gmsh mesh file PyFR/examples/couette_flow_2d/couette_flow_2d.msh into couette_flow_2d/

       4. Run pyfr to covert the Gmsh mesh file into a PyFR mesh file called couette_flow_2d.pyfrm:

             pyfr import couette_flow_2d.msh couette_flow_2d.pyfrm

       5. Run pyfr to solve the Navier-Stokes equations on the mesh, generating a series of PyFR solution  files
          called couette_flow_2d-*.pyfrs:

             pyfr run -b cuda -p couette_flow_2d.pyfrm couette_flow_2d.ini

       6. Run  pyfr  on  the solution file couette_flow_2d-040.pyfrs converting it into an unstructured VTK file
          called couette_flow_2d-040.vtu. Note that in order to visualise the high-order data,  each  high-order
          element  is  sub-divided  into smaller linear elements. The level of sub-division is controlled by the
          integer at the end of the command:

             pyfr export couette_flow_2d.pyfrm couette_flow_2d-040.pyfrs couette_flow_2d-040.vtu -d 4

       7. Visualise the unstructured VTK file in Paraview
         [image: couette flow] [image] Colour map of steady-state density distribution..UNINDENT

   Example --- 2D Euler Vortex
       Proceed with the following steps to run a parallel 2D Euler vortex simulation on a structured mesh:

       1. Create a working directory called euler_vortex_2d/

       2. Copy the configuration file PyFR/examples/euler_vortex_2d/euler_vortex_2d.ini into euler_vortex_2d/

       3. Copy the Gmsh file PyFR/examples/euler_vortex_2d/euler_vortex_2d.msh into euler_vortex_2d/

       4. Run pyfr to convert the Gmsh mesh file into a PyFR mesh file called euler_vortex_2d.pyfrm:

             pyfr import euler_vortex_2d.msh euler_vortex_2d.pyfrm

       5. Run pyfr to partition the PyFR mesh file into two pieces:

             pyfr partition 2 euler_vortex_2d.pyfrm .

       6. Run pyfr to solve the Euler equations on the mesh, generating a series of PyFR solution  files  called
          euler_vortex_2d*.pyfrs:

             mpirun -n 2 pyfr run -b cuda -p euler_vortex_2d.pyfrm euler_vortex_2d.ini

       7. Run  pyfr on the solution file euler_vortex_2d-100.0.pyfrs converting it into an unstructured VTK file
          called euler_vortex_2d-100.0.vtu. Note that in order to visualise the high-order data, each high-order
          element  is  sub-divided  into smaller linear elements. The level of sub-division is controlled by the
          integer at the end of the command:

             pyfr export euler_vortex_2d.pyfrm euler_vortex_2d-100.0.pyfrs euler_vortex_2d-100.0.vtu -d 4

       8. Visualise the unstructured VTK file in Paraview
         [image: euler vortex] [image] Colour map of density distribution at 100 time units..UNINDENT

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)

              _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)

              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)

              _abc_impl = <_abc_data object>

              _accept_step(dt, idxcurr, err=None)

              _add(*args)

              _controller_needs_errest

              _errest(x, y, z)

              _get_axnpby_kerns(n, subdims=None)

              _get_errest_kerns()

              _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)

              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)

              _abc_impl = <_abc_data object>

              _accept_step(dt, idxcurr)

              _add(*args)

              _dual_time_source()

              _get_axnpby_kerns(n, subdims=None)

              _get_errest_kerns()

              _get_gndofs()

              _get_kernels(name, nargs, **kwargs)

              _get_plugins()

              _get_reg_banks(nreg)

              _prepare_reg_banks(*bidxes)

              _resid(x, y)

              _source_regidx

              _stepper_nfevals

              _stepper_nregs

              _stepper_order

              _stepper_regidx

              advance_to(t)

              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:

   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)

              _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)

              stepper_name = 'euler'

       class pyfr.integrators.std.steppers.StdRK4Stepper(*args, **kwargs)

              _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)

              stepper_name = 'rk4'

       class pyfr.integrators.std.steppers.StdRK34Stepper(*args, **kwargs)

              _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)

              _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)

              _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)

              stepper_name = 'tvd-rk3'

       class pyfr.integrators.dual.steppers.DualBDF2Stepper(*args, **kwargs)

              _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)

              _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)

              _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)

              _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)

       class pyfr.integrators.dual.pseudosteppers.DualPseudoTVDRK3Stepper(*args, **kwargs)

              _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)

       class pyfr.integrators.dual.pseudosteppers.DualPseudoEulerStepper(*args, **kwargs)

              _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)

       Types of Stepper 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.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)

              _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)

              _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:

   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)

              _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)

              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)

              _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)

              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:

   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)

              _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)

              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)

              _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)

              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:

   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)

              _abc_impl = <_abc_data object>

              _malloc_impl(nbytes)

              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)

              _abc_impl = <_abc_data object>

              _malloc_impl(nbytes)

              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)

              _abc_impl = <_abc_data object>

              _malloc_impl(nbytes)

              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)

              _abc_impl = <_abc_data object>

              _malloc_impl(nbytes)

              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:

   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)

              _abc_impl = <_abc_data object>

              _build_arglst(dims, argn, argt, argdict)

              _build_kernel(name, src, argtypes)

              _instantiate_kernel(dims, fun, arglst)

              _render_kernel(name, mod, tplargs)

              kernel_generator_cls
                     alias of pyfr.backends.cuda.generator.CUDAKernelGenerator

              register(mod)

       class pyfr.backends.mic.provider.MICPointwiseKernelProvider(backend)

              _abc_impl = <_abc_data object>

              _build_arglst(dims, argn, argt, argdict)

              _build_kernel(name, src, argtypes, restype=None)

              _instantiate_kernel(dims, fun, arglst)

              _render_kernel(name, mod, tplargs)

              kernel_generator_cls
                     alias of pyfr.backends.mic.generator.MICKernelGenerator

              register(mod)

       class pyfr.backends.opencl.provider.OpenCLPointwiseKernelProvider(backend)

              _abc_impl = <_abc_data object>

              _build_arglst(dims, argn, argt, argdict)

              _build_kernel(name, src, argtypes)

              _instantiate_kernel(dims, fun, arglst)

              _render_kernel(name, mod, tplargs)

              kernel_generator_cls
                     alias of pyfr.backends.opencl.generator.OpenCLKernelGenerator

              register(mod)

       class pyfr.backends.openmp.provider.OpenMPPointwiseKernelProvider(backend)

              _abc_impl = <_abc_data object>

              _build_arglst(dims, argn, argt, argdict)

              _build_kernel(name, src, argtypes, restype=None)

              _instantiate_kernel(dims, fun, arglst)

              _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:

   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)

              _abc_impl = <_abc_data object>

              _deref_arg_array_1d(arg)

              _deref_arg_array_2d(arg)

              _deref_arg_view(arg)

              _render_body(body)

              _render_spec()

              argspec()

              needs_ldim(arg)

              render()

       class pyfr.backends.mic.generator.MICKernelGenerator(name, ndim, args, body, fpdtype)

              _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()

              argspec()

              needs_ldim(arg)

              render()

       class pyfr.backends.opencl.generator.OpenCLKernelGenerator(*args, **kwargs)

              _abc_impl = <_abc_data object>

              _deref_arg_array_1d(arg)

              _deref_arg_array_2d(arg)

              _deref_arg_view(arg)

              _render_body(body)

              _render_spec()

              argspec()

              needs_ldim(arg)

              render()

       class pyfr.backends.openmp.generator.OpenMPKernelGenerator(name, ndim, args, body, fpdtype)

              _abc_impl = <_abc_data object>

              _deref_arg_array_1d(arg)

              _deref_arg_array_2d(arg)

              _deref_arg_view(arg)

              _render_body(body)

              _render_spec()

              argspec()

              needs_ldim(arg)

              render()

       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

       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

INDICES AND TABLES

       • genindex

       • modindex

       • search

AUTHOR

       Imperial College London

COPYRIGHT

       2013-2019, Imperial College London