Generator Scripts

Writing generator scripts is the primary usage idiom of pystencils-sfg. A generator script is a Python script, say kernels.py, which contains pystencils-sfg code at the top level that, when executed, emits source code to a pair of files kernels.h and kernels.cpp. This guide describes how to write such a generator script, its structure, and how it can be used to generate code.

Anatomy

The code generation process in a generator script is controlled by the SourceFileGenerator context manager. It configures the code generator by combining configuration options from the environment (e.g. a CMake build system) with options specified in the script, and infers the names of the output files from the script’s name. It then returns a composer to the user, which provides a convenient interface for constructing the source files.

To start, place the following code in a Python script, e.g. kernels.py:

from pystencilssfg import SourceFileGenerator

with SourceFileGenerator() as sfg:
    pass

The source file is constructed within the context manager’s managed region. During execution of the script, when the region ends, a header/source file pair kernels.h and kernels.cpp will be written to disk next to your script. Execute the script as-is and inspect the generated files, which will of course still be empty:

Generated Files

#pragma once

#include <cstdint>

#define RESTRICT __restrict__


#include "kernels.h"



#define FUNC_PREFIX inline

/*************************************************************************************
 *                                Kernels
*************************************************************************************/



namespace kernels {


} // namespace kernels



/*************************************************************************************
 *                                Functions
*************************************************************************************/



/*************************************************************************************
 *                                Class Methods
*************************************************************************************/


Using the Composer

The object sfg constructed in above snippet is an instance of SfgComposer. The composer is the central part of the user front-end of pystencils-sfg. It provides an interface for constructing source files that closely mimics C++ syntactic structures within Python. Here is an overview of its various functions:

Includes and Definitions

With SfgComposer.include, the code generator can be instructed to include header files. As in C++, you can use the <> delimiters for system headers, and omit them for project headers.

from pystencilssfg import SourceFileGenerator

with SourceFileGenerator() as sfg:
    sfg.include("<vector>")
    sfg.include("<span>")
    sfg.include("custom_header.hpp")

#pragma once

#include "custom_header.hpp"
#include <span>
#include <vector>
#include <cstdint>

#define RESTRICT __restrict__


#include "kernels.h"



#define FUNC_PREFIX inline

/*************************************************************************************
 *                                Kernels
*************************************************************************************/



namespace kernels {


} // namespace kernels



/*************************************************************************************
 *                                Functions
*************************************************************************************/



/*************************************************************************************
 *                                Class Methods
*************************************************************************************/


Adding Kernels

pystencils-generated kernels are managed in kernel namespaces. The default kernel namespace is called kernels and is available via sfg.kernels. Adding an existing pystencils AST, or creating one from a list of assignments, is possible through kernels.add and kernels.create. The latter is a wrapper around pystencils.create_kernel. Both functions return a kernel handle through which the kernel can be accessed, e.g. for calling it in a function.

To access other kernel namespaces than the default one, the sfg.kernel_namespace method can be used.

from pystencilssfg import SourceFileGenerator

import pystencils as ps
import sympy as sp

with SourceFileGenerator() as sfg:
    #   Define a copy kernel
    src, dst = ps.fields("src, dst: [1D]")
    c = sp.Symbol("c")

    @ps.kernel
    def scale():
        dst.center @= c * src.center()

    #   Add it to the file
    scale_kernel = sfg.kernels.create(scale, "scale")

#pragma once

#include <cstdint>

#define RESTRICT __restrict__


#include "kernels.h"

#include <math.h>

#define FUNC_PREFIX inline

/*************************************************************************************
 *                                Kernels
*************************************************************************************/



namespace kernels {

FUNC_PREFIX void scale (double * const  _data_dst, double * const  _data_src, const int64_t _size_dst_0, const int64_t _stride_dst_0, const int64_t _stride_src_0, const double c)
{
   for(int64_t ctr_0 = 0LL; ctr_0 < _size_dst_0; ctr_0 += 1LL)
   {
      _data_dst[ctr_0 * _stride_dst_0] = c * _data_src[ctr_0 * _stride_src_0];
   }

}

} // namespace kernels



/*************************************************************************************
 *                                Functions
*************************************************************************************/



/*************************************************************************************
 *                                Class Methods
*************************************************************************************/


Building Functions

Through the composer, you can define free functions in your generated C++ file. These may contain arbitrary code; their primary intended task however is to wrap kernel calls with the necessary boilerplate code to integrate them into a framework. The composer provides an interface for constructing functions that tries to mimic the look of the generated C++ code. Use sfg.function to create a function, and sfg.call to call a kernel:

    #   ... see above ...
    sfg.function("scale_kernel")(
        sfg.call(scale_kernel)
    )

#pragma once

#include <cstdint>

#define RESTRICT __restrict__


 void scale_kernel ( double * const  _data_dst, double * const  _data_src, const int64_t _size_dst_0, const int64_t _stride_dst_0, const int64_t _stride_src_0, const double c );


#include "kernels.h"

#include <math.h>

#define FUNC_PREFIX inline

/*************************************************************************************
 *                                Kernels
*************************************************************************************/



namespace kernels {

FUNC_PREFIX void scale (double * const  _data_dst, double * const  _data_src, const int64_t _size_dst_0, const int64_t _stride_dst_0, const int64_t _stride_src_0, const double c)
{
   for(int64_t ctr_0 = 0LL; ctr_0 < _size_dst_0; ctr_0 += 1LL)
   {
      _data_dst[ctr_0 * _stride_dst_0] = c * _data_src[ctr_0 * _stride_src_0];
   }

}

} // namespace kernels



/*************************************************************************************
 *                                Functions
*************************************************************************************/



  void scale_kernel (double * const  _data_dst, double * const  _data_src, const int64_t _size_dst_0, const int64_t _stride_dst_0, const int64_t _stride_src_0, const double c){
  kernels::scale(_data_dst, _data_src, _size_dst_0, _stride_dst_0, _stride_src_0, c);}



/*************************************************************************************
 *                                Class Methods
*************************************************************************************/


Note the special syntax: To mimic the look of a C++ function, the composer uses a sequence of two calls to construct the function.

The function body can furthermore be populated with code to embedd the generated kernel into the target C++ application. If you examine the generated files of the previous example, you will notice that your function scale_kernel has lots of raw pointers and integer indices in its interface. We can wrap those up into proper C++ data structures, such as, for example, std::span or std::vector, like this:

    import pystencilssfg.lang.cpp.std as std

    sfg.include("<span>")
    
    sfg.function("scale_kernel")(
        sfg.map_field(src, std.vector(src)),
        sfg.map_field(dst, std.span(dst)),
        sfg.call(scale_kernel)
    )

#pragma once

#include <cstdint>
#include <vector>
#include <span>

#define RESTRICT __restrict__


 void scale_kernel ( const double c,  std::span< double > & dst,  std::vector< double > & src );


#include "kernels.h"

#include <math.h>

#define FUNC_PREFIX inline

/*************************************************************************************
 *                                Kernels
*************************************************************************************/



namespace kernels {

FUNC_PREFIX void scale (double * const  _data_dst, double * const  _data_src, const int64_t _size_dst_0, const int64_t _stride_dst_0, const int64_t _stride_src_0, const double c)
{
   for(int64_t ctr_0 = 0LL; ctr_0 < _size_dst_0; ctr_0 += 1LL)
   {
      _data_dst[ctr_0 * _stride_dst_0] = c * _data_src[ctr_0 * _stride_src_0];
   }

}

} // namespace kernels



/*************************************************************************************
 *                                Functions
*************************************************************************************/



  void scale_kernel (const double c,  std::span< double > & dst,  std::vector< double > & src){
  double * const  _data_src { src.data() };
  const int64_t _stride_src_0 { 1 };
  double * const  _data_dst { dst.data() };
  const int64_t _size_dst_0 { dst.size() };
  const int64_t _stride_dst_0 { 1 };
  kernels::scale(_data_dst, _data_src, _size_dst_0, _stride_dst_0, _stride_src_0, c);}



/*************************************************************************************
 *                                Class Methods
*************************************************************************************/


If you now inspect the generated code, you will see that the interface of your function is considerably simplified. Also, all the necessary code was added to its body to extract the low-level information required by the actual kernel from the data structures.

The sfg.map_field API can be used to map pystencils fields to a variety of different data structures. The pystencils-sfg provides modelling support for a number of C++ standard library classes (see pystencilssfg.lang.cpp.std). It also provides the necessary infrastructure for modelling the data structures of any C++ framework in a similar manner.