psyclone.transformations

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.

Classes

class psyclone.transformations.ACCEnterDataTrans

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

Inheritance

Inheritance diagram of ACCEnterDataTrans
apply(sched, options=None)

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)

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.ACCDataTrans

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)

Inheritance

Inheritance diagram of ACCDataTrans
apply(node, options=None)

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)

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

Parameters:
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.ACCKernelsTrans

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)

Inheritance

Inheritance diagram of ACCKernelsTrans
apply(node, options=None)

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)

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

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

Inheritance

Inheritance diagram of ACCLoopTrans
apply(node, options=None)

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)

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

Inheritance

Inheritance diagram of ACCParallelTrans
apply(target_nodes, options=None)

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)

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

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)

Inheritance

Inheritance diagram of ACCRoutineTrans
apply(node, options=None)

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:
property name
Returns:

the name of this transformation class.

Return type:

str

validate(node, options=None)

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

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

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

  • options["force"] (bool) – whether to allow routines with CodeBlocks to run on the GPU.

Raises:
class psyclone.transformations.ColourTrans

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

Inheritance

Inheritance diagram of ColourTrans
apply(node, options=None)

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

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

Inheritance

Inheritance diagram of Dynamo0p3AsyncHaloExchangeTrans
apply(node, options=None)

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)

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

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

Inheritance

Inheritance diagram of Dynamo0p3ColourTrans
apply(node, options=None)

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

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

Inheritance

Inheritance diagram of Dynamo0p3KernelConstTrans
apply(node, options=None)

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)

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

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’.

Inheritance

Inheritance diagram of Dynamo0p3OMPLoopTrans
apply(node, options=None)

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)

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

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.

Inheritance

Inheritance diagram of Dynamo0p3RedundantComputationTrans
apply(loop, options=None)

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)

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:
class psyclone.transformations.DynamoOMPParallelLoopTrans(omp_directive='do', omp_schedule='static')

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’.

Inheritance

Inheritance diagram of DynamoOMPParallelLoopTrans
validate(node, options=None)

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:
class psyclone.transformations.GOceanOMPLoopTrans(omp_directive='do', omp_schedule='static')

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’.

Inheritance

Inheritance diagram of GOceanOMPLoopTrans
validate(node, options=None)

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

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’.

Inheritance

Inheritance diagram of GOceanOMPParallelLoopTrans
apply(node, options=None)

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

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.

Inheritance

Inheritance diagram of KernelImportsToArguments
apply(node, options=None)

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)

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

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.

Inheritance

Inheritance diagram of MoveTrans
apply(node, location, options=None)

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:
property name

Returns the name of this transformation as a string.

validate(node, location, options=None)

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:
class psyclone.transformations.OMPLoopTrans(omp_directive='do', omp_schedule='auto')

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

Inheritance

Inheritance diagram of OMPLoopTrans
apply(node, options=None)

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

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

Inheritance

Inheritance diagram of OMPMasterTrans
property name
Returns:

the name of this transformation as a string.

Return type:

str

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

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

Inheritance

Inheritance diagram of OMPParallelLoopTrans
apply(node, options=None)

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

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

Inheritance

Inheritance diagram of OMPParallelTrans
property name
Returns:

the name of this transformation as a string.

Return type:

str

validate(node_list, options=None)

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)

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

Inheritance

Inheritance diagram of OMPSingleTrans
apply(node_list, options=None)

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

Base class for transformations that create a parallel region.

Inheritance

Inheritance diagram of ParallelRegionTrans
apply(target_nodes, options=None)

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)

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