API

The generator module

This module provides the PSyclone ‘main’ routine which is intended to be driven from the bin/psyclone executable script. ‘main’ takes an algorithm file as input and produces modified algorithm code and generated PSy code. A function, ‘generate’, is also provided which has the same functionality as ‘main’ but can be called from within another Python program.

psyclone.generator.generate(filename, api='', kernel_paths=None, script_name=None, line_length=False, distributed_memory=None, kern_out_path='', kern_naming='multiple')[source]

Takes a PSyclone algorithm specification as input and outputs the associated generated algorithm and psy codes suitable for compiling with the specified kernel(s) and support infrastructure. Uses the parse.algorithm.parse() function to parse the algorithm specification, the psyGen.PSy class to generate the PSy code and the alg_gen.Alg class to generate the modified algorithm code.

Parameters:
  • filename (str) – the file containing the algorithm specification.

  • api (str) – the name of the API to use. Defaults to empty string.

  • kernel_paths (Optional[List[str]]) – the directories from which to recursively search for the files containing the kernel source (if different from the location of the algorithm specification). Defaults to None.

  • script_name (str) – a script file that can apply optimisations to the PSy layer (can be a path to a file or a filename that relies on the PYTHONPATH to find the module). Defaults to None.

  • line_length (bool) – a logical flag specifying whether we care about line lengths being longer than 132 characters. If so, the input (algorithm and kernel) code is checked to make sure that it conforms. The default is False.

  • distributed_memory (bool) – a logical flag specifying whether to generate distributed memory code. The default is set in the ‘config.py’ file.

  • kern_out_path (str) – directory to which to write transformed kernel code. Defaults to empty string.

  • kern_naming (bool) – the scheme to use when re-naming transformed kernels. Defaults to “multiple”.

Returns:

2-tuple containing the fparser1 AST for the algorithm code and the fparser1 AST or a string (for NEMO) of the psy code.

Return type:

Tuple[fparser.one.block_statements.BeginSource, fparser.one.block_statements.Module] | Tuple[fparser.one.block_statements.BeginSource, str]

Raises:
  • GenerationError – if an invalid API is specified.

  • GenerationError – if an invalid kernel-renaming scheme is specified.

  • GenerationError – if there is an error raising the PSyIR to domain-specific PSyIR.

  • GenerationError – if a kernel functor is not named in a use statement.

  • IOError – if the filename or search path do not exist.

  • NoInvokesError – if no invokes are found in the algorithm file.

For example:

>>> from psyclone.generator import generate
>>> alg, psy = generate("algspec.f90")
>>> alg, psy = generate("algspec.f90", kernel_paths=["src/kernels"])
>>> alg, psy = generate("algspec.f90", script_name="optimise.py")
>>> alg, psy = generate("algspec.f90", line_length=True)
>>> alg, psy = generate("algspec.f90", distributed_memory=False)

The parse module

Module that uses the Fortran parser fparser2 to parse PSyclone-conformant Algorithm code.

psyclone.parse.algorithm.parse(alg_filename, api='', invoke_name='invoke', kernel_paths=None, line_length=False)[source]

Takes a PSyclone conformant algorithm file as input and outputs a parse tree of the code contained therein and an object containing information about the ‘invoke’ calls in the algorithm file and any associated kernels within the invoke calls.

Parameters:
  • alg_filename (str) – the file containing the algorithm specification.

  • api (str) – the PSyclone API to use when parsing the code. Defaults to empty string.

  • invoke_name (str) – the expected name of the invocation calls in the algorithm code. Defaults to “invoke”.

  • kernel_paths (list of str or NoneType) – the paths to search for kernel source files (if different from the location of the algorithm source). Defaults to None.

  • line_length (bool) – a logical flag specifying whether we care about line lengths being longer than 132 characters. If so, the input (algorithm and kernel) code is checked to make sure that it conforms and an error raised if not. The default is False.

Returns:

2-tuple consisting of the fparser2 parse tree of the Algorithm file and an object holding details of the invokes found.

Return type:

(fparser.two.Fortran2003.Program, psyclone.parse.FileInfo)

For example:

>>> from psyclone.parse.algorithm import parse
>>> ast, info = parse(SOURCE_FILE)

The transformations module

This module provides the various transformations that can be applied to PSyIR nodes. There are both general and API-specific transformation classes in this module where the latter typically apply API-specific checks before calling the base class for the actual transformation.

class psyclone.transformations.ACCDataTrans[source]

Add an OpenACC data region around a list of nodes in the PSyIR. COPYIN, COPYOUT and COPY clauses are added as required.

For example:

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.psyGen import PSyFactory
>>> api = "nemo"
>>> ast, invokeInfo = parse(NEMO_SOURCE_FILE, api=api)
>>> psy = PSyFactory(api).create(invokeInfo)
>>>
>>> from psyclone.transformations import ACCKernelsTrans, ACCDataTrans
>>> ktrans = ACCKernelsTrans()
>>> dtrans = ACCDataTrans()
>>>
>>> schedule = psy.invokes.get('tra_adv').schedule
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
>>>
>>> # Add a kernels construct for execution on the device
>>> kernels = schedule.children[9]
>>> ktrans.apply(kernels)
>>>
>>> # Enclose the kernels in a data construct
>>> kernels = schedule.children[9]
>>> dtrans.apply(kernels)
apply(node, options=None)[source]

Put the supplied node or list of nodes within an OpenACC data region.

Parameters:
  • node ((list of) psyclone.psyir.nodes.Node) – the PSyIR node(s) to enclose in the data region.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

property name
Returns:

the name of this transformation.

Return type:

str

validate(nodes, options)[source]

Check that we can safely add a data region around the supplied list of nodes.

Parameters:
  • nodes (List[psyclone.psyir.nodes.Node] | psyclone.psyir.nodes.Node) – the proposed node(s) to enclose in a data region.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

Raises:
  • TransformationError – if the Schedule to which the nodes belong already has an ‘enter data’ directive.

  • TransformationError – if any of the nodes are themselves data directives.

  • TransformationError – if an array of structures needs to be deep copied (this is not currently supported).

class psyclone.transformations.ACCEnterDataTrans[source]

Adds an OpenACC “enter data” directive to a Schedule. For example:

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.psyGen import PSyFactory
>>> api = "gocean1.0"
>>> ast, invokeInfo = parse(GOCEAN_SOURCE_FILE, api=api)
>>> psy = PSyFactory(api).create(invokeInfo)
>>>
>>> from psyclone.transformations import         ACCEnterDataTrans, ACCLoopTrans, ACCParallelTrans
>>> dtrans = ACCEnterDataTrans()
>>> ltrans = ACCLoopTrans()
>>> ptrans = ACCParallelTrans()
>>>
>>> schedule = psy.invokes.get('invoke_0').schedule
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
>>>
>>> # Apply the OpenACC Loop transformation to *every* loop in the schedule
>>> for child in schedule.children[:]:
...     ltrans.apply(child)
>>>
>>> # Enclose all of these loops within a single OpenACC parallel region
>>> ptrans.apply(schedule)
>>>
>>> # Add an enter data directive
>>> dtrans.apply(schedule)
>>>
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
apply(sched, options=None)[source]

Adds an OpenACC “enter data” directive to the invoke associated with the supplied Schedule. Any fields accessed by OpenACC kernels within this schedule will be added to this data region in order to ensure they remain on the target device.

Parameters:
  • sched (sub-class of psyclone.psyir.nodes.Schedule) – schedule to which to add an “enter data” directive.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

property name
Returns:

the name of this transformation.

Return type:

str

validate(sched, options=None)[source]

Check that we can safely apply the OpenACC enter-data transformation to the supplied Schedule.

Parameters:
  • sched (sub-class of psyclone.psyir.nodes.Schedule) – Schedule to which to add an “enter data” directive.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

Raises:

TransformationError – if passed something that is not a (subclass of) psyclone.psyir.nodes.Schedule.

class psyclone.transformations.ACCKernelsTrans[source]

Enclose a sub-set of nodes from a Schedule within an OpenACC kernels region (i.e. within “!$acc kernels” … “!$acc end kernels” directives).

For example:

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.psyGen import PSyFactory
>>> api = "nemo"
>>> ast, invokeInfo = parse(NEMO_SOURCE_FILE, api=api)
>>> psy = PSyFactory(api).create(invokeInfo)
>>>
>>> from psyclone.transformations import ACCKernelsTrans
>>> ktrans = ACCKernelsTrans()
>>>
>>> schedule = psy.invokes.get('tra_adv').schedule
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
>>> kernels = schedule.children[9]
>>> # Transform the kernel
>>> ktrans.apply(kernels)
apply(node, options=None)[source]

Enclose the supplied list of PSyIR nodes within an OpenACC Kernels region.

Parameters:
  • node ((a list of) psyclone.psyir.nodes.Node) – a node or list of nodes in the PSyIR to enclose.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

  • options["default_present"] (bool) – whether or not the kernels region should have the ‘default present’ attribute (indicating that data is already on the accelerator). When using managed memory this option should be False.

property name
Returns:

the name of this transformation class.

Return type:

str

validate(nodes, options)[source]

Check that we can safely enclose the supplied node or list of nodes within OpenACC kernels … end kernels directives.

Parameters:
  • nodes ((list of) psyclone.psyir.nodes.Node) – the proposed PSyIR node or nodes to enclose in the kernels region.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

Raises:
  • NotImplementedError – if the supplied Nodes belong to a GOInvokeSchedule.

  • TransformationError – if there are no Loops within the proposed region.

class psyclone.transformations.ACCLoopTrans[source]

Adds an OpenACC loop directive to a loop. This directive must be within the scope of some OpenACC Parallel region (at code-generation time).

For example:

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.parse.utils import ParseError
>>> from psyclone.psyGen import PSyFactory
>>> from psyclone.errors import GenerationError
>>> api = "gocean1.0"
>>> ast, invokeInfo = parse(GOCEAN_SOURCE_FILE, api=api)
>>> psy = PSyFactory(api).create(invokeInfo)
>>>
>>> from psyclone.psyGen import TransInfo
>>> t = TransInfo()
>>> ltrans = t.get_trans_name('ACCLoopTrans')
>>> rtrans = t.get_trans_name('ACCParallelTrans')
>>>
>>> schedule = psy.invokes.get('invoke_0').schedule
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
>>>
>>> # Apply the OpenACC Loop transformation to *every* loop in the schedule
>>> for child in schedule.children[:]:
...     ltrans.apply(child)
>>>
>>> # Enclose all of these loops within a single OpenACC parallel region
>>> rtrans.apply(schedule)
>>>
apply(node, options=None)[source]

Apply the ACCLoop transformation to the specified node. This node must be a Loop since this transformation corresponds to inserting a directive immediately before a loop, e.g.:

!$ACC LOOP
do ...
   ...
end do

At code-generation time (when psyclone.psyir.nodes.ACCLoopDirective.gen_code() is called), this node must be within (i.e. a child of) a PARALLEL region.

Parameters:
  • node (psyclone.psyir.nodes.Loop) – the supplied node to which we will apply the Loop transformation.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

  • options["collapse"] (int) – number of nested loops to collapse.

  • options["independent"] (bool) – whether to add the “independent” clause to the directive (not strictly necessary within PARALLEL regions).

  • options["sequential"] (bool) – whether to add the “seq” clause to the directive.

  • options["gang"] (bool) – whether to add the “gang” clause to the directive.

  • options["vector"] (bool) – whether to add the “vector” clause to the directive.

class psyclone.transformations.ACCParallelTrans(default_present=True)[source]

Create an OpenACC parallel region by inserting an ‘acc parallel’ directive.

>>> from psyclone.psyGen import TransInfo
>>> from psyclone.psyir.frontend.fortran import FortranReader
>>> from psyclone.psyir.backend.fortran import FortranWriter
>>> from psyclone.psyir.nodes import Loop
>>> psyir = FortranReader().psyir_from_source("""
... program do_loop
...     real, dimension(10) :: A
...     integer i
...     do i = 1, 10
...       A(i) = i
...     end do
... end program do_loop
... """)
>>> ptrans = TransInfo().get_trans_name('ACCParallelTrans')
>>>
>>> # Enclose the loop within a OpenACC PARALLEL region
>>> ptrans.apply(psyir.walk(Loop))
>>> print(FortranWriter()(psyir))
program do_loop
  real, dimension(10) :: a
  integer :: i

  !$acc parallel default(present)
  do i = 1, 10, 1
    a(i) = i
  enddo
  !$acc end parallel

end program do_loop
apply(target_nodes, options=None)[source]

Encapsulate given nodes with the ACCParallelDirective.

Parameters:
  • target_nodes (psyclone.psyir.nodes.Node | List[psyclone.psyir.nodes.Node]) – a single Node or a list of Nodes.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

  • options["node-type-check"] (bool) – this flag controls if the type of the nodes enclosed in the region should be tested to avoid using unsupported nodes inside a region.

  • options["default_present"] (bool) – this flag controls if the inserted directive should include the default_present clause.

validate(node_list, options=None)[source]

Validate this transformation.

Parameters:
  • node_list (psyclone.psyir.nodes.Node | List[psyclone.psyir.nodes.Node]) – a single Node or a list of Nodes.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

  • options["node-type-check"] (bool) – this flag controls if the type of the nodes enclosed in the region should be tested to avoid using unsupported nodes inside a region.

  • options["default_present"] (bool) – this flag controls if the inserted directive should include the default_present clause.

class psyclone.transformations.ACCRoutineTrans[source]

Transform a kernel or routine by adding a “!$acc routine” directive (causing it to be compiled for the OpenACC accelerator device). For example:

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.psyGen import PSyFactory
>>> api = "gocean1.0"
>>> ast, invokeInfo = parse(GOCEAN_SOURCE_FILE, api=api)
>>> psy = PSyFactory(api).create(invokeInfo)
>>>
>>> from psyclone.transformations import ACCRoutineTrans
>>> rtrans = ACCRoutineTrans()
>>>
>>> schedule = psy.invokes.get('invoke_0').schedule
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
>>> kern = schedule.children[0].children[0].children[0]
>>> # Transform the kernel
>>> rtrans.apply(kern)
apply(node, options=None)[source]

Add the ‘!$acc routine’ OpenACC directive into the code of the supplied Kernel (in a PSyKAl API such as GOcean or LFRic) or directly in the supplied Routine.

Parameters:
  • node (psyclone.psyGen.Kern or psyclone.psyir.nodes.Routine) – the kernel call or routine implementation to transform.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

property name
Returns:

the name of this transformation class.

Return type:

str

validate(node, options=None)[source]

Perform checks that the supplied kernel or routine can be transformed.

Parameters:
  • node (psyclone.psyGen.Kern | psyclone.psyir.nodes.Routine) – the kernel which is the target of the transformation.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

Raises:
  • TransformationError – if the node is not a kernel or a routine.

  • TransformationError – if the target is a built-in kernel.

  • TransformationError – if it is a kernel but without an associated PSyIR.

  • TransformationError – if any of the symbols in the kernel are accessed via a module use statement.

  • TransformationError – if the kernel contains any calls to other routines.

class psyclone.transformations.ColourTrans[source]

Apply a colouring transformation to a loop (in order to permit a subsequent parallelisation over colours). For example:

>>> invoke = ...
>>> schedule = invoke.schedule
>>>
>>> ctrans = ColourTrans()
>>>
>>> # Colour all of the loops
>>> for child in schedule.children:
>>>     ctrans.apply(child)
>>>
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
apply(node, options=None)[source]

Converts the Loop represented by node into a nested loop where the outer loop is over colours and the inner loop is over cells of that colour.

Parameters:
  • node (psyclone.psyir.nodes.Loop) – the loop to transform.

  • options (Optional[Dict[str, Any]]) – options for the transformation.

class psyclone.transformations.Dynamo0p3AsyncHaloExchangeTrans[source]

Splits a synchronous halo exchange into a halo exchange start and halo exchange end. For example:

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.psyGen import PSyFactory
>>> api = "dynamo0.3"
>>> ast, invokeInfo = parse("file.f90", api=api)
>>> psy=PSyFactory(api).create(invokeInfo)
>>> schedule = psy.invokes.get('invoke_0').schedule
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
>>>
>>> from psyclone.transformations import Dynamo0p3AsyncHaloExchangeTrans
>>> trans = Dynamo0p3AsyncHaloExchangeTrans()
>>> trans.apply(schedule.children[0])
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
apply(node, options=None)[source]

Transforms a synchronous halo exchange, represented by a HaloExchange node, into an asynchronous halo exchange, represented by HaloExchangeStart and HaloExchangeEnd nodes.

Parameters:
  • node (psyclone.psygen.HaloExchange) – a synchronous haloexchange node.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

property name
Returns:

the name of this transformation as a string.

Return type:

str

validate(node, options)[source]

Internal method to check whether the node is valid for this transformation.

Parameters:
  • node (psyclone.psygen.HaloExchange) – a synchronous Halo Exchange node

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

Raises:

TransformationError – if the node argument is not a HaloExchange (or subclass thereof)

class psyclone.transformations.Dynamo0p3ColourTrans[source]

Split a Dynamo 0.3 loop over cells into colours so that it can be parallelised. For example:

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.psyGen import PSyFactory
>>> import transformations
>>> import os
>>> import pytest
>>>
>>> TEST_API = "dynamo0.3"
>>> _,info=parse(os.path.join(os.path.dirname(os.path.abspath(__file__)),
>>>              "tests", "test_files", "dynamo0p3",
>>>              "4.6_multikernel_invokes.f90"),
>>>              api=TEST_API)
>>> psy = PSyFactory(TEST_API).create(info)
>>> invoke = psy.invokes.get('invoke_0')
>>> schedule = invoke.schedule
>>>
>>> ctrans = Dynamo0p3ColourTrans()
>>> otrans = DynamoOMPParallelLoopTrans()
>>>
>>> # Colour all of the loops
>>> for child in schedule.children:
>>>     ctrans.apply(child)
>>>
>>> # Then apply OpenMP to each of the colour loops
>>> for child in schedule.children:
>>>     otrans.apply(child.children[0])
>>>
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())

Colouring in the LFRic (Dynamo 0.3) API is subject to the following rules:

  • Only kernels which operate on ‘CELL_COLUMN’s and which increment a field on a continuous function space require colouring. Kernels that update a field on a discontinuous function space will cause this transformation to raise an exception. Kernels that only write to a field on a continuous function space also do not require colouring but are permitted.

  • A kernel may have at most one field with ‘GH_INC’ access.

  • A separate colour map will be required for each field that is coloured (if an invoke contains >1 kernel call).

apply(node, options=None)[source]

Performs Dynamo0.3-specific error checking and then uses the parent class to convert the Loop represented by node into a nested loop where the outer loop is over colours and the inner loop is over cells of that colour.

Parameters:
  • node (psyclone.domain.lfric.LFRicLoop) – the loop to transform.

  • options – a dictionary with options for transformations. :type options: Optional[Dict[str, Any]]

class psyclone.transformations.Dynamo0p3KernelConstTrans[source]

Modifies a kernel so that the number of dofs, number of layers and number of quadrature points are fixed in the kernel rather than being passed in by argument.

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.psyGen import PSyFactory
>>> api = "dynamo0.3"
>>> ast, invokeInfo = parse("file.f90", api=api)
>>> psy=PSyFactory(api).create(invokeInfo)
>>> schedule = psy.invokes.get('invoke_0').schedule
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
>>>
>>> from psyclone.transformations import Dynamo0p3KernelConstTrans
>>> trans = Dynamo0p3KernelConstTrans()
>>> for kernel in schedule.coded_kernels():
>>>     trans.apply(kernel, number_of_layers=150)
>>>     kernel_schedule = kernel.get_kernel_schedule()
>>>     # Uncomment the following line to see a text view of the
>>>     # symbol table
>>>     # print(kernel_schedule.symbol_table.view())
apply(node, options=None)[source]

Transforms a kernel so that the values for the number of degrees of freedom (if a valid value for the element_order arg is provided), the number of quadrature points (if the quadrature arg is set to True) and the number of layers (if a valid value for the number_of_layers arg is provided) are constant in a kernel rather than being passed in by argument.

The “cellshape”, “element_order” and “number_of_layers” arguments are provided to mirror the namelist values that are input into an LFRic model when it is run.

Quadrature support is currently limited to XYoZ in ths transformation. In the case of XYoZ the number of quadrature points (for horizontal and vertical) are set to the element_order + 3 in the LFRic infrastructure so their value is derived.

Parameters:
  • node (psyclone.domain.lfric.LFRicKern) – a kernel node.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

  • options["cellshape"] (str) – the shape of the cells. This is provided as it helps determine the number of dofs a field has for a particular function space. Currently only “quadrilateral” is supported which is also the default value.

  • options["element_order"] (int) – the order of the cell. In combination with cellshape, this determines the number of dofs a field has for a particular function space. If it is set to None (the default) then the dofs values are not set as constants in the kernel, otherwise they are.

  • options["number_of_layers"] (int) – the number of vertical layers in the LFRic model mesh used for this particular run. If this is set to None (the default) then the nlayers value is not set as a constant in the kernel, otherwise it is.

  • options["quadrature"] (bool) – whether the number of quadrature points values are set as constants in the kernel (True) or not (False). The default is False.

property name
Returns:

the name of this transformation as a string.

Return type:

str

validate(node, options=None)[source]

This method checks whether the input arguments are valid for this transformation.

Parameters:
  • node (psyclone.domain.lfric.LFRicKern) – a dynamo 0.3 kernel node.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

  • options["cellshape"] (str) – the shape of the elements/cells.

  • options["element_order"] (int) – the order of the elements/cells.

  • options["number_of_layers"] (int) – the number of layers to use.

  • options["quadrature"] (bool) – whether quadrature dimension sizes should or shouldn’t be set as constants in a kernel.

Raises:

TransformationError – if the node argument is not a dynamo 0.3 kernel, the cellshape argument is not set to “quadrilateral”, the element_order argument is not a 0 or a positive integer, the number of layers argument is not a positive integer, the quadrature argument is not a boolean, neither element order nor number of layers arguments are set (as the transformation would then do nothing), or the quadrature argument is True but the element order is not provided (as the former needs the latter).

class psyclone.transformations.Dynamo0p3OMPLoopTrans(omp_schedule='static')[source]

LFRic (Dynamo 0.3) specific orphan OpenMP loop transformation. Adds Dynamo-specific validity checks.

Parameters:

omp_schedule (str) – the OpenMP schedule to use. Must be one of ‘runtime’, ‘static’, ‘dynamic’, ‘guided’ or ‘auto’. Defaults to ‘static’.

apply(node, options=None)[source]

Apply LFRic (Dynamo 0.3) specific OMPLoopTrans.

Parameters:
  • node (psyclone.psyir.nodes.Node) – the Node in the Schedule to check.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations and validation.

  • options["reprod"] (bool) – indicating whether reproducible reductions should be used. By default the value from the config file will be used.

validate(node, options=None)[source]

Perform LFRic (Dynamo 0.3) specific loop validity checks for the OMPLoopTrans.

Parameters:
  • node (psyclone.psyir.nodes.Node) – the Node in the Schedule to check

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations and validation.

  • options["reprod"] (bool) – indicating whether reproducible reductions should be used. By default the value from the config file will be used.

Raises:

TransformationError – if an OMP loop transform would create incorrect code.

class psyclone.transformations.Dynamo0p3RedundantComputationTrans[source]

This transformation allows the user to modify a loop’s bounds so that redundant computation will be performed. Redundant computation can result in halo exchanges being modified, new halo exchanges being added or existing halo exchanges being removed.

  • This transformation should be performed before any parallelisation transformations (e.g. for OpenMP) to the loop in question and will raise an exception if this is not the case.

  • This transformation can not be applied to a loop containing a reduction and will again raise an exception if this is the case.

  • This transformation can only be used to add redundant computation to a loop, not to remove it.

  • This transformation allows a loop that is already performing redundant computation to be modified, but only if the depth is increased.

apply(loop, options=None)[source]

Apply the redundant computation transformation to the loop loop. This transformation can be applied to loops iterating over ‘cells or ‘dofs’. if depth is set to a value then the value will be the depth of the field’s halo over which redundant computation will be performed. If depth is not set to a value then redundant computation will be performed to the full depth of the field’s halo.

Parameters:
  • loop (psyclone.psyGen.LFRicLoop) – the loop that we are transforming.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

  • options["depth"] (int) – the depth of the stencil. Defaults to None.

validate(node, options=None)[source]

Perform various checks to ensure that it is valid to apply the RedundantComputation transformation to the supplied node

Parameters:
  • node (psyclone.psyir.nodes.Node) – the supplied node on which we are performing validity checks

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

  • options["depth"] (int) – the depth of the stencil if the value is provided and None if not.

Raises:
  • TransformationError – if the parent of the loop is a psyclone.psyir.nodes.Directive.

  • TransformationError – if the parent of the loop is not a psyclone.psyir.nodes.Loop or a psyclone.psyGen.DynInvokeSchedule.

  • TransformationError – if the parent of the loop is a psyclone.psyir.nodes.Loop but the original loop does not iterate over ‘colour’.

  • TransformationError – if the parent of the loop is a psyclone.psyir.nodes.Loop but the parent does not iterate over ‘colours’.

  • TransformationError – if the parent of the loop is a psyclone.psyir.nodes.Loop but the parent’s parent is not a psyclone.psyGen.DynInvokeSchedule.

  • TransformationError – if this transformation is applied when distributed memory is not switched on.

  • TransformationError – if the loop does not iterate over cells, dofs or colour.

  • TransformationError – if the transformation is setting the loop to the maximum halo depth but the loop already computes to the maximum halo depth.

  • TransformationError – if the transformation is setting the loop to the maximum halo depth but the loop contains a stencil access (as this would result in the field being accessed beyond the halo depth).

  • TransformationError – if the supplied depth value is not an integer.

  • TransformationError – if the supplied depth value is less than 1.

  • TransformationError – if the supplied depth value is not greater than 1 when a continuous loop is modified as this is the minimum valid value.

  • TransformationError – if the supplied depth value is not greater than the existing depth value, as we should not need to undo existing transformations.

  • TransformationError – if a depth value has been supplied but the loop has already been set to the maximum halo depth.

class psyclone.transformations.DynamoOMPParallelLoopTrans(omp_directive='do', omp_schedule='static')[source]

Dynamo-specific OpenMP loop transformation. Adds Dynamo specific validity checks. Actual transformation is done by the base class.

Parameters:
  • omp_directive (str) – choose which OpenMP loop directive to use. Defaults to “do”.

  • omp_schedule (str) – the OpenMP schedule to use. Must be one of ‘runtime’, ‘static’, ‘dynamic’, ‘guided’ or ‘auto’. Defaults to ‘static’.

validate(node, options=None)[source]

Perform LFRic-specific loop validity checks then call the validate method of the base class.

Parameters:
  • node (psyclone.psyir.nodes.Node) – the Node in the Schedule to check

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

Raises:
  • TransformationError – if the supplied Node is not a LFRicLoop.

  • TransformationError – if the associated loop requires colouring.

class psyclone.transformations.GOceanOMPLoopTrans(omp_directive='do', omp_schedule='static')[source]

GOcean-specific orphan OpenMP loop transformation. Adds GOcean specific validity checks (that the node is either an inner or outer Loop).

Parameters:
  • omp_directive (str) – choose which OpenMP loop directive to use. Defaults to “do”.

  • omp_schedule (str) – the OpenMP schedule to use. Must be one of ‘runtime’, ‘static’, ‘dynamic’, ‘guided’ or ‘auto’. Defaults to ‘static’.

validate(node, options=None)[source]

Checks that the supplied node is a valid target for parallelisation using OMP directives.

Parameters:
  • node (psyclone.psyir.nodes.Loop) – the candidate loop for parallelising using OMP Do.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

Raises:

TransformationError – if the loop_type of the supplied Loop is not “inner” or “outer”.

class psyclone.transformations.GOceanOMPParallelLoopTrans(omp_directive='do', omp_schedule='static')[source]

GOcean specific OpenMP Do loop transformation. Adds GOcean specific validity checks (that supplied Loop is an inner or outer loop). Actual transformation is done by base class.

param str omp_directive:

choose which OpenMP loop directive to use. Defaults to “do”.

param str omp_schedule:

the OpenMP schedule to use. Must be one of ‘runtime’, ‘static’, ‘dynamic’, ‘guided’ or ‘auto’. Defaults to ‘static’.

apply(node, options=None)[source]

Perform GOcean-specific loop validity checks then call OMPParallelLoopTrans.apply().

Parameters:
  • node (psyclone.psyir.nodes.Loop) – a Loop node from an AST.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations and validation.

Raises:

TransformationError – if the supplied node is not an inner or outer loop.

class psyclone.transformations.KernelImportsToArguments[source]

Transformation that removes any accesses of imported data from the supplied kernel and places them in the caller. The values/references are then passed by argument into the kernel.

apply(node, options=None)[source]

Convert the imported variables used inside the kernel into arguments and modify the InvokeSchedule to pass the same imported variables to the kernel call.

Parameters:
  • node (psyclone.psyGen.CodedKern) – a kernel call.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

property name
Returns:

the name of this transformation.

Return type:

str

validate(node, options=None)[source]

Check that the supplied node is a valid target for this transformation.

Parameters:
  • node (psyclone.psyGen.CodedKern) – the PSyIR node to validate.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

Raises:
  • TransformationError – if the supplied node is not a CodedKern.

  • TransformationError – if this transformation is not applied to a Gocean API Invoke.

  • TransformationError – if the supplied kernel contains wildcard imports of symbols from one or more containers (e.g. a USE without an ONLY clause in Fortran).

class psyclone.transformations.MoveTrans[source]

Provides a transformation to move a node in the tree. For example:

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.psyGen import PSyFactory
>>> ast,invokeInfo=parse("dynamo.F90")
>>> psy=PSyFactory("dynamo0.3").create(invokeInfo)
>>> schedule=psy.invokes.get('invoke_v3_kernel_type').schedule
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
>>>
>>> from psyclone.transformations import MoveTrans
>>> trans=MoveTrans()
>>> trans.apply(schedule.children[0], schedule.children[2],
...             options = {"position":"after")
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())

Nodes may only be moved to a new location with the same parent and must not break any dependencies otherwise an exception is raised.

apply(node, location, options=None)[source]

Move the node represented by node before location location (which is also a node) by default and after if the optional position argument is set to ‘after’.

Parameters:
  • node (psyclone.psyir.nodes.Node) – the node to be moved.

  • location (psyclone.psyir.nodes.Node) – node before or after which the given node should be moved.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

  • options["position"] (str) – either ‘before’ or ‘after’.

Raises:
  • TransformationError – if the given node is not an instance of psyclone.psyir.nodes.Node

  • TransformationError – if the location is not valid.

property name

Returns the name of this transformation as a string.

validate(node, location, options=None)[source]

validity checks for input arguments.

Parameters:
  • node (psyclone.psyir.nodes.Node) – the node to be moved.

  • location (psyclone.psyir.nodes.Node) – node before or after which the given node should be moved.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

  • options["position"] (str) – either ‘before’ or ‘after’.

Raises:
  • TransformationError – if the given node is not an instance of psyclone.psyir.nodes.Node

  • TransformationError – if the location is not valid.

class psyclone.transformations.OMPLoopTrans(omp_directive='do', omp_schedule='auto')[source]

Adds an OpenMP directive to parallelise this loop. It can insert different directives such as “omp do/for”, “omp parallel do/for”, “omp teams distribute parallel do/for” or “omp loop” depending on the provided parameters. The OpenMP schedule to use can also be specified, but this will be ignored in case of the “omp loop” (as the ‘schedule’ clause is not valid for this specific directive). The configuration-defined ‘reprod’ parameter also specifies whether a manual reproducible reproduction is to be used. Note, reproducible in this case means obtaining the same results with the same number of OpenMP threads, not for different numbers of OpenMP threads.

Parameters:
  • omp_schedule (str) – the OpenMP schedule to use. Defaults to ‘auto’.

  • omp_directive (str) – choose which OpenMP loop directive to use. Defaults to “omp do”

For example:

>>> from psyclone.psyir.frontend.fortran import FortranReader
>>> from psyclone.psyir.backend.fortran import FortranWriter
>>> from psyclone.psyir.nodes import Loop
>>> from psyclone.transformations import OMPLoopTrans, OMPParallelTrans
>>>
>>> psyir = FortranReader().psyir_from_source("""
...     subroutine my_subroutine()
...         integer, dimension(10, 10) :: A
...         integer :: i
...         integer :: j
...         do i = 1, 10
...             do j = 1, 10
...                 A(i, j) = 0
...             end do
...         end do
...     end subroutine
...     """)
>>> loop = psyir.walk(Loop)[0]
>>> omplooptrans1 = OMPLoopTrans(omp_schedule="dynamic",
...                              omp_directive="paralleldo")
>>> omplooptrans1.apply(loop)
>>> print(FortranWriter()(psyir))
subroutine my_subroutine()
  integer, dimension(10,10) :: a
  integer :: i
  integer :: j

  !$omp parallel do default(shared), private(i,j), schedule(dynamic)
  do i = 1, 10, 1
    do j = 1, 10, 1
      a(i,j) = 0
    enddo
  enddo
  !$omp end parallel do

end subroutine my_subroutine
apply(node, options=None)[source]

Apply the OMPLoopTrans transformation to the specified PSyIR Loop.

Parameters:
  • node (psyclone.psyir.nodes.Node) – the supplied node to which we will apply the OMPLoopTrans transformation

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations and validation.

  • options["reprod"] (bool) – indicating whether reproducible reductions should be used. By default the value from the config file will be used.

property omp_directive
Returns:

the type of OMP directive that this transformation will insert.

Return type:

str

property omp_schedule
Returns:

the OpenMP schedule that will be specified by this transformation.

Return type:

str

class psyclone.transformations.OMPMasterTrans[source]

Create an OpenMP MASTER region by inserting directives. The most likely use case for this transformation is to wrap around task-based transformations. Note that adding this directive requires a parent OpenMP parallel region (which can be inserted by OMPParallelTrans), otherwise it will produce an error in generation-time.

For example:

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.psyGen import PSyFactory
>>> api = "gocean1.0"
>>> ast, invokeInfo = parse(GOCEAN_SOURCE_FILE, api=api)
>>> psy = PSyFactory(api).create(invokeInfo)
>>>
>>> from psyclone.transformations import OMPParallelTrans, OMPMasterTrans
>>> mastertrans = OMPMasterTrans()
>>> paralleltrans = OMPParallelTrans()
>>>
>>> schedule = psy.invokes.get('invoke_0').schedule
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
>>>
>>> # Enclose all of these loops within a single OpenMP
>>> # MASTER region
>>> mastertrans.apply(schedule.children)
>>> # Enclose all of these loops within a single OpenMP
>>> # PARALLEL region
>>> paralleltrans.apply(schedule.children)
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
property name
Returns:

the name of this transformation as a string.

Return type:

str

class psyclone.transformations.OMPParallelLoopTrans(omp_directive='do', omp_schedule='auto')[source]

Adds an OpenMP PARALLEL DO directive to a loop.

For example:

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.psyGen import PSyFactory
>>> ast, invokeInfo = parse("dynamo.F90")
>>> psy = PSyFactory("dynamo0.3").create(invokeInfo)
>>> schedule = psy.invokes.get('invoke_v3_kernel_type').schedule
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
>>>
>>> from psyclone.transformations import OMPParallelLoopTrans
>>> trans = OMPParallelLoopTrans()
>>> trans.apply(schedule.children[0])
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
apply(node, options=None)[source]

Apply an OMPParallelLoop Transformation to the supplied node (which must be a Loop). In the generated code this corresponds to wrapping the Loop with directives:

!$OMP PARALLEL DO ...
do ...
  ...
end do
!$OMP END PARALLEL DO
Parameters:
  • node (psyclone.f2pygen.DoGen) – the node (loop) to which to apply the transformation.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations and validation.

class psyclone.transformations.OMPParallelTrans[source]

Create an OpenMP PARALLEL region by inserting directives. For example:

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.parse.utils import ParseError
>>> from psyclone.psyGen import PSyFactory
>>> from psyclone.errors import GenerationError
>>> api = "gocean1.0"
>>> ast, invokeInfo = parse(GOCEAN_SOURCE_FILE, api=api)
>>> psy = PSyFactory(api).create(invokeInfo)
>>>
>>> from psyclone.psyGen import TransInfo
>>> t = TransInfo()
>>> ltrans = t.get_trans_name('GOceanOMPLoopTrans')
>>> rtrans = t.get_trans_name('OMPParallelTrans')
>>>
>>> schedule = psy.invokes.get('invoke_0').schedule
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
>>>
>>> # Apply the OpenMP Loop transformation to *every* loop
>>> # in the schedule
>>> for child in schedule.children:
>>>     ltrans.apply(child)
>>>
>>> # Enclose all of these loops within a single OpenMP
>>> # PARALLEL region
>>> rtrans.apply(schedule.children)
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
property name
Returns:

the name of this transformation as a string.

Return type:

str

validate(node_list, options=None)[source]

Perform OpenMP-specific validation checks.

Parameters:
  • node_list (list of psyclone.psyir.nodes.Node) – list of Nodes to put within parallel region.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

  • options["node-type-check"] (bool) – this flag controls if the type of the nodes enclosed in the region should be tested to avoid using unsupported nodes inside a region.

Raises:

TransformationError – if the target Nodes are already within some OMP parallel region.

class psyclone.transformations.OMPSingleTrans(nowait=False)[source]

Create an OpenMP SINGLE region by inserting directives. The most likely use case for this transformation is to wrap around task-based transformations. The parent region for this should usually also be a OMPParallelTrans.

Parameters:

nowait (bool) – whether to apply a nowait clause to this transformation. The default value is False

For example:

>>> from psyclone.parse.algorithm import parse
>>> from psyclone.psyGen import PSyFactory
>>> api = "gocean1.0"
>>> ast, invokeInfo = parse(GOCEAN_SOURCE_FILE, api=api)
>>> psy = PSyFactory(api).create(invokeInfo)
>>>
>>> from psyclone.transformations import OMPParallelTrans, OMPSingleTrans
>>> singletrans = OMPSingleTrans()
>>> paralleltrans = OMPParallelTrans()
>>>
>>> schedule = psy.invokes.get('invoke_0').schedule
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
>>>
>>> # Enclose all of these loops within a single OpenMP
>>> # SINGLE region
>>> singletrans.apply(schedule.children)
>>> # Enclose all of these loops within a single OpenMP
>>> # PARALLEL region
>>> paralleltrans.apply(schedule.children)
>>> # Uncomment the following line to see a text view of the schedule
>>> # print(schedule.view())
apply(node_list, options=None)[source]

Apply the OMPSingleTrans transformation to the specified node in a Schedule.

At code-generation time this node must be within (i.e. a child of) an OpenMP PARALLEL region. Code generation happens when OMPLoopDirective.gen_code() is called, or when the PSyIR tree is given to a backend.

If the keyword “nowait” is specified in the options, it will cause a nowait clause to be added if it is set to True, otherwise no clause will be added.

Parameters:
  • node_list ((a list of) psyclone.psyir.nodes.Node) – the supplied node or node list to which we will apply the OMPSingleTrans transformation

  • options (Optional[Dict[str, Any]]) – a list with options for transformations and validation.

  • options["nowait"] (bool) – indicating whether or not to use a nowait clause on this single region.

property name
Returns:

the name of this transformation.

Return type:

str

property omp_nowait
Returns:

whether or not this Single region uses a nowait clause to remove the end barrier.

Return type:

bool

class psyclone.transformations.ParallelRegionTrans[source]

Base class for transformations that create a parallel region.

apply(target_nodes, options=None)[source]

Apply this transformation to a subset of the nodes within a schedule - i.e. enclose the specified Loops in the schedule within a single parallel region.

Parameters:
  • target_nodes ((list of) psyclone.psyir.nodes.Node) – a single Node or a list of Nodes.

  • options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.

  • options["node-type-check"] (bool) – this flag controls if the type of the nodes enclosed in the region should be tested to avoid using unsupported nodes inside a region.

validate(node_list, options=None)[source]

Check that the supplied list of Nodes are eligible to be put inside a parallel region.

Parameters:
  • node_list (list) – list of nodes to put into a parallel region

  • options – a dictionary with options for transformations. :type options: Optional[Dict[str, Any]]

  • options["node-type-check"] (bool) – this flag controls whether or not the type of the nodes enclosed in the region should be tested to avoid using unsupported nodes inside a region.

Raises:
  • TransformationError – if the supplied node is an InvokeSchedule rather than being within an InvokeSchedule.

  • TransformationError – if the supplied nodes are not all children of the same parent (siblings).

The psyGen module

This module provides generic support for PSyclone’s PSy code optimisation and generation. The classes in this method need to be specialised for a particular API and implementation.

class psyclone.psyGen.PSy(invoke_info)[source]

Base class to help manage and generate PSy code for a single algorithm file. Takes the invocation information output from the function parse.algorithm.parse() as its input and stores this in a way suitable for optimisation and code generation.

Parameters:

invoke_info (psyclone.parse.algorithm.FileInfo) – An object containing the required invocation information for code optimisation and generation. Produced by the function parse.algorithm.parse().

For example:

>>> from psyclone.parse.algorithm import parse
>>> ast, info = parse("argspec.F90")
>>> from psyclone.psyGen import PSyFactory
>>> api = "..."
>>> psy = PSyFactory(api).create(info)
>>> print(psy.gen)
property container
Returns:

the container associated with this PSy object

Return type:

psyclone.psyir.nodes.Container

abstract property gen

Abstract base class for code generation function.

Returns:

root node of generated Fortran AST.

Return type:

psyclone.psyir.nodes.Node

property invokes
Returns:

the list of invokes.

Return type:

psyclone.psyGen.Invokes or derived class

property name
Returns:

the name of the PSy object.

Return type:

str

The alg_gen module

This module provides the Alg class and supporting exception-handling to translate the original algorithm file into one that can be compiled and linked with the generated PSy code.

class psyclone.alg_gen.Alg(parse_tree, psy, invoke_name='invoke')[source]

Generate a modified algorithm code for a single algorithm specification. Takes the parse tree of the algorithm specification output from the function psyclone.parse.algorithm.parse() and an instance of the psyGen.PSy class as input. The latter allows consistent names to be generated between the algorithm (callng) and psy (callee) layers.

For example:

>>> from psyclone.algorithm.parse import parse
>>> parse_tree, info = parse("argspec.F90")
>>> from psyclone.psyGen import PSy
>>> psy = PSy(info)
>>> from psyclone.alg_gen import Alg
>>> alg = Alg(parse_tree, psy)
>>> print(alg.gen)
Parameters:
  • parse_tree (fparser.two.utils.Base) – an object containing a parse tree of the algorithm specification which was produced by the function psyclone.parse.algorithm.parse(). Assumes the algorithm will be parsed by fparser2 and expects a valid program unit, program, module, subroutine or function.

  • psy (psyclone.psyGen.PSy) – an object containing information about the PSy layer.

  • invoke_name (str) – the name that the algorithm layer uses to indicate an invoke call. This is an optional argument that defaults to the name “invoke”.

property gen

Modifies and returns the algorithm code. ‘invoke’ calls are replaced with calls to the corresponding PSy-layer routines and the USE statements for the kernels that were referenced by each ‘invoke’ are removed.

Returns:

the modified algorithm specification as an fparser2 parse tree.

Return type:

fparser.two.utils.Base

Raises:

NoInvokesError – if no ‘invoke()’ calls are found.

The line_length module

Provides support for breaking long fortran lines into smaller ones to allow the code to conform to the maximum line length limits (132 for f90 free format is the default)

class psyclone.line_length.FortLineLength(line_length=132)[source]

This class take a free format fortran code as a string and line wraps any lines that are larger than the specified line length

property length

returns the maximum allowed line length

long_lines(fortran_in)[source]

returns true if at least one of the lines in the input code is longer than the allowed length. Otherwise returns false

process(fortran_in)[source]

Processes unlimited line-length Fortran code into Fortran code with long lines wrapped appropriately.

Parameters:

fortran_in (str) – Fortran code to be line wrapped.

Returns:

line wrapped Fortran code.

Return type:

str