psyclone.psyir.transformations#
Transformation module, containing all generic (API independent) transformations and base classes.
Submodules#
psyclone.psyir.transformations.acc_kernels_transpsyclone.psyir.transformations.acc_update_transpsyclone.psyir.transformations.allarrayaccess2loop_transpsyclone.psyir.transformations.arrayaccess2loop_transpsyclone.psyir.transformations.arrayassignment2loops_transpsyclone.psyir.transformations.async_trans_mixinpsyclone.psyir.transformations.chunk_loop_transpsyclone.psyir.transformations.debug_checksum_transpsyclone.psyir.transformations.extract_transpsyclone.psyir.transformations.fold_conditional_return_expressions_transpsyclone.psyir.transformations.hoist_local_arrays_transpsyclone.psyir.transformations.hoist_loop_bound_expr_transpsyclone.psyir.transformations.hoist_transpsyclone.psyir.transformations.increase_rank_loop_arrays_transpsyclone.psyir.transformations.inline_transpsyclone.psyir.transformations.intrinsics- Submodules
psyclone.psyir.transformations.intrinsics.abs2code_transpsyclone.psyir.transformations.intrinsics.array_reduction_base_transpsyclone.psyir.transformations.intrinsics.dotproduct2code_transpsyclone.psyir.transformations.intrinsics.intrinsic2code_transpsyclone.psyir.transformations.intrinsics.matmul2code_transpsyclone.psyir.transformations.intrinsics.max2code_transpsyclone.psyir.transformations.intrinsics.maxval2loop_transpsyclone.psyir.transformations.intrinsics.min2code_transpsyclone.psyir.transformations.intrinsics.minormax2code_transpsyclone.psyir.transformations.intrinsics.minval2loop_transpsyclone.psyir.transformations.intrinsics.product2loop_transpsyclone.psyir.transformations.intrinsics.sign2code_transpsyclone.psyir.transformations.intrinsics.sum2loop_trans
- Submodules
psyclone.psyir.transformations.loop_fuse_transpsyclone.psyir.transformations.loop_swap_transpsyclone.psyir.transformations.loop_tiling_2d_transpsyclone.psyir.transformations.loop_tiling_transpsyclone.psyir.transformations.loop_transpsyclone.psyir.transformations.mark_routine_for_gpu_mixinpsyclone.psyir.transformations.omp_declare_target_transpsyclone.psyir.transformations.omp_loop_transpsyclone.psyir.transformations.omp_minimise_sync_transpsyclone.psyir.transformations.omp_target_transpsyclone.psyir.transformations.omp_task_transpsyclone.psyir.transformations.omp_taskloop_transpsyclone.psyir.transformations.omp_taskwait_transpsyclone.psyir.transformations.parallel_loop_transpsyclone.psyir.transformations.parallel_region_transpsyclone.psyir.transformations.profile_transpsyclone.psyir.transformations.psy_data_transpsyclone.psyir.transformations.read_only_verify_transpsyclone.psyir.transformations.reference2arrayrange_transpsyclone.psyir.transformations.region_transpsyclone.psyir.transformations.replace_induction_variables_transpsyclone.psyir.transformations.replace_reference_by_literal_transpsyclone.psyir.transformations.scalarisation_transpsyclone.psyir.transformations.transformation_errorpsyclone.psyir.transformations.value_range_check_trans
Classes#
ACCKernelsTrans: Enclose a sub-set of nodes from a Schedule within an OpenACC kernelsACCUpdateTrans: Examines the supplied Schedule and adds OpenACC “update” directivesAllArrayAccess2LoopTrans: Provides a transformation from a PSyIR Assignment containingArrayAccess2LoopTrans: Provides a transformation to transform a constant index access toArrayAssignment2LoopsTrans: Provides a transformation from a PSyIR Array Range to a PSyIRDebugChecksumTrans: Creates a set of checksums (written via print) for all written arraysChunkLoopTrans: Apply a chunking transformation to a loop (in order to permit aExtractTrans: This transformation inserts an ExtractNode or a node derivedFoldConditionalReturnExpressionsTrans: Provides a transformation that folds conditional expressions with onlyHoistLocalArraysTrans: This transformation takes a Routine and promotes any local arrays toHoistLoopBoundExprTrans: This transformation moves complex bounds expressions out of the loopHoistTrans: This transformation takes an assignment and moves it outside ofIncreaseRankLoopArraysTrans: This transformation takes a loop and a list of arrays accessed insideInlineTrans: This transformation takes a Call (which may have a return value)Abs2CodeTrans: Provides a transformation from a PSyIR ABS Operator node toDotProduct2CodeTrans: Provides a transformation from a PSyIR DOT_PRODUCT Operator node toMatmul2CodeTrans: Provides a transformation from a PSyIR MATMUL Operator node toMax2CodeTrans: Provides a transformation from a PSyIR MAX Intrinsic node toMin2CodeTrans: Provides a transformation from a PSyIR MIN Intrinsic node toSign2CodeTrans: Provides a transformation from a PSyIR SIGN intrinsic node toSum2LoopTrans: Provides a transformation from a PSyIR SUM IntrinsicCall node to anLoopFuseTrans: Provides a generic loop-fuse transformation to two Nodes in theLoopSwapTrans: Provides a loop-swap transformation, e.g.:LoopTilingTrans: Apply a loop tiling transformation to a loop. For example:LoopTiling2DTrans: Apply a 2D loop tiling transformation to a loop. This is a specialLoopTrans: This abstract class is a base class for all transformations that actMaxval2LoopTrans: Provides a transformation from a PSyIR MAXVAL IntrinsicCall node toMinval2LoopTrans: Provides a transformation from a PSyIR MINVAL IntrinsicCall node toOMPLoopTrans: Adds an OpenMP directive to parallelise this loop. It can insert differentOMPMinimiseSyncTrans: Attempts to remove OMPTaskwaitDirective orOMPTargetTrans: Adds an OpenMP target directive to a region of code.OMPTaskTrans: Apply an OpenMP Task Transformation to a Loop. The Loop mustOMPTaskwaitTrans: Adds zero or more OpenMP Taskwait directives to an OMP parallel region.ParallelLoopTrans: Adds an abstract directive (it needs to be specified by sub-classing thisProduct2LoopTrans: Provides a transformation from a PSyIR PRODUCT IntrinsicCall node toProfileTrans: Create a profile region around a list of statements. ForPSyDataTrans: Create a PSyData region around a list of statements. ForReadOnlyVerifyTrans: This transformation inserts a ReadOnlyVerifyNode or a node derivedReference2ArrayRangeTrans: Provides a transformation from PSyIR Array Notation (a reference toRegionTrans: This abstract class is a base class for all transformations that actReplaceInductionVariablesTrans: Move all supported induction variables out of the loop, and replaceReplaceReferenceByLiteralTrans: This transformation takes a psyir Routine and replace all Reference psyirValueRangeCheckTrans: This transformation inserts a ValueRangeCheckNode into the PSyIR of aParallelRegionTrans: Base class for transformations that create a parallel region.OMPTaskloopTrans: Adds an OpenMP taskloop directive to a loop. Only one of grainsize orOMPDeclareTargetTrans: Adds an OpenMP declare target directive to the specified routine.
- class psyclone.psyir.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.psyir.frontend import FortranReader >>> psyir = FortranReader().psyir_from_source(NEMO_SOURCE_FILE) >>> >>> from psyclone.psyir.transformations import ACCKernelsTrans >>> ktrans = ACCKernelsTrans() >>> >>> schedule = psyir.children[0] >>> # 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

- apply(node, options={})[source]#
Enclose the supplied list of PSyIR nodes within an OpenACC Kernels region.
- Parameters:
node (
Union[Node,List[Node]]) – a node or list of nodes in the PSyIR to enclose.options (
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.
options["disable_loop_check"] (bool) – whether to disable the check that the supplied region contains 1 or more loops. Default is False (i.e. the check is enabled).
options["async_queue"] (Union[bool,
psyclone.psyir.nodes.DataNode]) – whether or not to add the ‘async’ clause to the new directive and if so, which queue to associate it with. True to enable for the default queue or a queue value specified with an int or PSyIR expression.options["allow_string"] (bool) – whether to allow the transformation on assignments involving character types. Defaults to False.
options["verbose"] (bool) – log the reason the validation failed, at the moment with a comment in the provided PSyIR node.
- static check_async_queue(nodes, async_queue)[source]#
Common function to check that all parent data directives have the same async queue.
- Parameters:
node – the nodes in the PSyIR to enclose.
async_queue (
Union[bool,int,Reference]) – The async queue to expect in ancestors.
- Raises:
TypeError – if the supplied queue is of the wrong type.
TransformationError – if the supplied queue does not match that specified by any ancestor nodes.
- 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 (
Union[Node,List[Node]]) – the proposed PSyIR node or nodes to enclose in the kernels region.options (
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.
options["disable_loop_check"] (bool) – whether to disable the check that the supplied region contains 1 or more loops. Default is False (i.e. the check is enabled).
options["async_queue"] (Union[bool,
psyclone.psyir.nodes.DataNode]) – whether or not to add the ‘async’ clause to the new directive and if so, which queue to associate it with. True to enable for the default queue or a queue value specified with an int or PSyIR expression.options["allow_string"] (bool) – whether to allow the transformation on assignments involving character types. Defaults to False.
options["verbose"] (bool) – log the reason the validation failed, at the moment with a comment in the provided PSyIR node.
- Raises:
NotImplementedError – if the supplied Nodes belong to a GOInvokeSchedule.
TransformationError – if there is an access to an assumed-size character variable within the region.
TransformationError – if the proposed region contains a call to a routine that is not available on the accelerator.
TransformationError – if there are no Loops within the proposed region and options[“disable_loop_check”] is not True.
TransformationError – if any assignments in the region contain a character type child and options[“allow_string”] is not True.
- Return type:
None
- class psyclone.psyir.transformations.ACCUpdateTrans[source]#
Examines the supplied Schedule and adds OpenACC “update” directives for any data accessed outside of a kernels or parallel region. For example:
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Routine >>> from psyclone.psyir.transformations import ACCUpdateTrans >>> >>> code = """ ... subroutine run_it() ... real :: a(10) ... a(:) = 7.0 ... end subroutine run_it""" >>> >>> psyir = FortranReader().psyir_from_source(code) >>> routine = psyir.walk(Routine)[0] >>> >>> # Add update directives >>> uptrans = ACCUpdateTrans() >>> uptrans.apply(routine) >>> >>> # Uncomment to see a text view of the new routine schedule >>> # print(routine.view()) >>> print(FortranWriter()(routine)) subroutine run_it() real, dimension(10) :: a !$acc update if_present host(a) a(:) = 7.0 !$acc update if_present device(a) end subroutine run_it
Inheritance

- apply(node, options=None)[source]#
Applies this transformation to the supplied Schedule. Identifies any regions of code outside of OpenACC regions and adds the necessary OpenACC update directives to ensure that the host and device copies of any variables are kept up-to-date.
- Parameters:
node (
psyclone.psyir.nodes.Schedule) – the Schedule that is to be transformed.options (Optional[Dict[str, Any]]) – any options to this transformation.
- validate(node, options=None)[source]#
Checks that it is valid to apply this transformation to the supplied schedule.
- Parameters:
node (
psyclone.psyir.nodes.Schedule) – the Schedule that is to be transformed.options (Optional[Dict[str, Any]]) – any options to this transformation.
- Raises:
TransformationError – if the supplied node is not a Schedule.
TransformationError – if the supplied node is within an OpenACC region.
- class psyclone.psyir.transformations.AllArrayAccess2LoopTrans[source]#
Provides a transformation from a PSyIR Assignment containing constant index accesses to an array into single-trip loops. For example:
>>> from psyclone.psyir.transformations import AllArrayAccess2LoopTrans >>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Assignment >>> code = ("program example\n" ... " real a(10,10), b(10,10)\n" ... " integer :: n\n" ... " a(1,n-1) = b(1,n-1)\n" ... "end program example\n") >>> psyir = FortranReader().psyir_from_source(code) >>> assignment = psyir.walk(Assignment)[0] >>> AllArrayAccess2LoopTrans().apply(assignment) >>> print(FortranWriter()(psyir)) program example real, dimension(10,10) :: a real, dimension(10,10) :: b integer :: n integer :: idx integer :: idx_1 do idx = 1, 1, 1 do idx_1 = n - 1, n - 1, 1 a(idx,idx_1) = b(idx,idx_1) enddo enddo end program example
Inheritance

- apply(node, options=None)[source]#
Apply the AllArrayAccess2Loop transformation if the supplied node is an Assignment with an Array Reference on its left-hand-side. Each constant array index access (i.e. one not containing a loop iterator or a range) is then transformed into an iterator and the assignment placed within a single-trip loop, subject to any constraints in the ArrayAccess2Loop transformation.
If any of the AllArrayAccess2Loop constraints are not satisfied for a loop index then this transformation does nothing for that index and silently moves to the next.
- Parameters:
node (
psyclone.psyir.nodes.Assignment) – an assignment.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations. This is an optional argument that defaults to None.
- property name#
- Returns:
the name of the transformation as a string.
- Return type:
str
- validate(node, options=None)[source]#
Perform any checks to ensure that it is valid to apply the AllArrayAccess2LoopTrans transformation to the supplied PSyIR Node.
- Parameters:
node (
psyclone.psyir.nodes.Node) – the node that is being checked.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations. This is an optional argument that defaults to None.
- class psyclone.psyir.transformations.ArrayAccess2LoopTrans[source]#
Provides a transformation to transform a constant index access to an array (i.e. one that does not contain a loop iterator) to a single trip loop. For example:
>>> from psyclone.psyir.transformations import ArrayAccess2LoopTrans >>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Assignment >>> code = ("program example\n" ... " real a(10)\n" ... " a(1) = 0.0\n" ... "end program example\n") >>> psyir = FortranReader().psyir_from_source(code) >>> assignment = psyir.walk(Assignment)[0] >>> ArrayAccess2LoopTrans().apply(assignment.lhs.children[0]) >>> print(FortranWriter()(psyir)) program example real, dimension(10) :: a integer :: ji do ji = 1, 1, 1 a(ji) = 0.0 enddo end program example
Inheritance

- apply(node, options=None)[source]#
Apply the ArrayAccess2Loop transformation if the supplied node is an access to an array index within an Array Reference that is on the left-hand-side of an Assignment node. The access must be a scalar (i.e. not a range) and must not include a loop variable (as we are transforming a single access to a loop).
If the constraints are satisfied then the array access is replaced with a loop iterator and a single trip loop. The new loop will be placed immediately around the assignment.
- Parameters:
node (
psyclone.psyir.nodes.Node) – an array index.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations. This is an optional argument that defaults to None.
- validate(node, options=None)[source]#
Perform various checks to ensure that it is valid to apply the ArrayAccess2LoopTrans transformation to the supplied PSyIR Node.
- Parameters:
node (
psyclone.psyir.nodes.Node) – the node that is being checked.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations. This is an optional argument that defaults to None.
- class psyclone.psyir.transformations.ArrayAssignment2LoopsTrans[source]#
Provides a transformation from a PSyIR Array Range to a PSyIR Loop. For example:
>>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Assignment >>> from psyclone.psyir.transformations import ArrayAssignment2LoopsTrans >>> code = """ ... subroutine sub() ... real :: tmp(10) ... tmp(:) = tmp(:) + 3 ... end subroutine sub""" >>> psyir = FortranReader().psyir_from_source(code) >>> assignment = psyir.walk(Assignment)[0] >>> trans = ArrayAssignment2LoopsTrans() >>> trans.apply(assignment) >>> print(psyir.debug_string()) subroutine sub() real, dimension(10) :: tmp integer :: idx do idx = LBOUND(tmp, dim=1), UBOUND(tmp, dim=1), 1 tmp(idx) = tmp(idx) + 3 enddo end subroutine sub
By default the transformation will reject character arrays, though this can be overriden by setting the ‘allow_string’ option to True. Note that PSyclone expresses syntax such as character(LEN=100) as UnsupportedFortranType, and this transformation will convert unknown or unsupported types to loops.
Inheritance

- apply(node, options=None)[source]#
Apply the transformation to the specified array assignment node. Each range node within the assignment is replaced with an explicit loop. The bounds of the loop are determined from the bounds of the array range on the left hand side of the assignment.
- Parameters:
node (
psyclone.psyir.nodes.Assignment) – an Assignment node.options["allow_string"] (bool) – whether to allow the transformation on a character type array range. Defaults to False.
options["verbose"] (bool) – log the reason the validation failed, at the moment with a comment in the provided PSyIR node.
- validate(node, options=None)[source]#
Perform various checks to ensure that it is valid to apply the ArrayAssignment2LoopsTrans transformation to the supplied PSyIR Node.
By default the validate function will throw an TransofmrationError on character arrays, though this can be overriden by setting the allow_string option to True. Note that PSyclone expresses syntax such as character(LEN=100) as UnsupportedFortranType, and this transformation will convert unknown or unsupported types to loops.
- Parameters:
node (
psyclone.psyir.nodes.Assignment) – the node that is being checked.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations
options["allow_string"] (bool) – whether to allow the transformation on a character type array range. Defaults to False.
options["verbose"] (bool) – log the reason the validation failed, at the moment with a comment in the provided PSyIR node.
- Raises:
TransformationError – if the node argument is not an Assignment.
TransformationError – if the node argument is an Assignment whose left hand side is not an ArrayReference.
TransformationError – if the node argument is an Assignment whose left hand side is an ArrayReference that does not have Range specifying the access to at least one of its dimensions.
TransformationError – if two or more of the loop ranges in the assignment are different or are not known to be the same.
TransformationError – if node contains a character type child and the allow_strings option is not set.
- static validate_no_char(node, message, options)[source]#
Check that there is no character variable accessed in the sub-tree with the supplied node at its root.
- Parameters:
node (
Node) – the root node to check for character assignments.message (
str) – the message to use if a character assignment is found.options (
dict) – any options that apply to this check.options["verbose"] (bool) – log the reason the validation failed, at the moment with a comment in the provided PSyIR node.
- Raises:
TransformationError – if the supplied node contains a child of character type.
- Return type:
None
- class psyclone.psyir.transformations.DebugChecksumTrans[source]#
Creates a set of checksums (written via print) for all written arrays inside the provided region.
Warning
This transformation generates code that will only work with .F90 suffixed files. When using this transformation make sure the output file declared to PSyclone is .F90 and not .f90.
For example:
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.transformations import DebugChecksumTrans
>>> psyir = FortranReader().psyir_from_source(""" ... subroutine mysubroutine() ... integer, dimension(10,10) :: A ... integer :: i ... integer :: j ... do i = 1, 10 ... do j = 1, 10 ... A(i,j) = A(i,k) + i-j ... end do ... end do ... end subroutine ... """) ... loop = psyir.children[0].children[0] ... DebugChecksumTrans().apply(loop) ... print(FortranWriter()(psyir)) subroutine mysubroutine() integer, dimension(10,10) :: a integer :: i integer :: j integer :: PSYCLONE_INTERNAL_line_ do i = 1, 10, 1 do j = 1, 10, 1 a(i,j) = a(i,j) + i - j enddo enddo PSYCLONE_INTERNAL_line_ = __LINE__ PRINT *, "PSyclone checksums from mysubroutine at line:", PSYCLONE_INTERNAL_line_ + 1 PRINT *, "a checksum", SUM(a) end subroutine mysubroutine
Inheritance

- class psyclone.psyir.transformations.ChunkLoopTrans[source]#
Apply a chunking transformation to a loop (in order to permit a chunked parallelisation). For example:
>>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Loop >>> from psyclone.psyir.transformations import ChunkLoopTrans >>> psyir = FortranReader().psyir_from_source(""" ... subroutine sub() ... integer :: ji, tmp(100) ... do ji=1, 100 ... tmp(ji) = 2 * ji ... enddo ... end subroutine sub""") >>> loop = psyir.walk(Loop)[0] >>> ChunkLoopTrans().apply(loop)
will generate:
subroutine sub() integer :: ji integer, dimension(100) :: tmp integer :: ji_el_inner integer :: ji_out_var do ji_out_var = 1, 100, 32 ji_el_inner = MIN(ji_out_var + (32 - 1), 100) do ji = ji_out_var, ji_el_inner, 1 tmp(ji) = 2 * ji enddo enddo end subroutine sub
Inheritance

- apply(node, options=None)[source]#
Converts the given Loop node into a nested loop where the outer loop is over chunks and the inner loop is over each individual element of the chunk.
- Parameters:
node (
psyclone.psyir.nodes.Loop) – the loop to transform.options (Optional[Dict[str, Any]]) – a dict with options for transformations.
options["chunksize"] (int) – The size to chunk over for this transformation. If not specified, the value 32 is used.
- validate(node, options=None)[source]#
Validates that the given Loop node can have a ChunkLoopTrans applied.
- Parameters:
node (
psyclone.psyir.nodes.Loop) – the loop to validate.options (Optional[Dict[str, Any]]) – a dict with options for transformation.
options["chunksize"] (int) – The size to chunk over for this transformation. If not specified, the value 32 is used.
- Raises:
TransformationError – if the supplied Loop has a step size which is not a constant value.
TransformationError – if the supplied Loop has a non-integer step size.
TransformationError – if the supplied Loop has a step size larger than the chosen chunk size.
TransformationError – if the supplied Loop is a chunked loop.
TransformationError – if the supplied Loop has a step size of 0.
TransformationError – if the supplied Loop writes to Loop variables inside the Loop body.
TransformationError – if the supplied Loop contains a CodeBlock node.
TransformationError – if an unsupported option has been provided.
TransformationError – if the provided tilesize is not a positive integer.
- class psyclone.psyir.transformations.ExtractTrans(node_class=<class 'psyclone.psyir.nodes.extract_node.ExtractNode'>)[source]#
This transformation inserts an ExtractNode or a node derived from ExtractNode into the PSyIR of a schedule. At code creation time this node will use the PSyData API to create code that can write the input and output parameters to a file. The node might also create a stand-alone driver program that can read the created file and then execute the instrumented region. Examples are given in the derived classes LFRicExtractTrans and GOceanExtractTrans.
After applying the transformation the Nodes marked for extraction are children of the ExtractNode. Nodes to extract can be individual constructs within an Invoke (e.g. Loops containing a Kernel or BuiltIn call) or entire Invokes. This functionality does not support distributed memory.
- Parameters:
node_class (
psyclone.psyir.nodes.ExtractNodeor derived class) – The Node class of which an instance will be inserted into the tree (defaults to ExtractNode), but can be any derived class.
Inheritance

- validate(node_list, options=None)[source]#
Performs validation checks specific to extract-based transformations.
- Parameters:
node_list (list of
psyclone.psyir.nodes.Node) – the list of Node(s) we are checking.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- Raises:
TransformationError – if transformation is applied to a Kernel or a BuiltIn call without its parent Loop.
TransformationError – if transformation is applied to a Loop without its parent Directive when optimisations are applied.
TransformationError – if transformation is applied to an orphaned Directive without its parent Directive.
- class psyclone.psyir.transformations.FoldConditionalReturnExpressionsTrans[source]#
Provides a transformation that folds conditional expressions with only a return statement inside so that the Return statement is moved to the end of the Routine and therefore it can be safely removed. This simplifies the control flow of the kernel to facilitate other transformations like kernel fusions. For example, the following code:
subroutine test(i) if (i < 5) then return endif if (i > 10) then return endif ! CODE end subroutine
will be transformed to:
subroutine test(i) if (.not.(i < 5)) then if (.not.(i > 10)) then ! CODE endif endif end subroutine
Inheritance

- apply(node, options=None)[source]#
Apply this transformation to the supplied node.
- Parameters:
node (
psyclone.psyir.nodes.Routine) – the node to transform.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- property name#
Returns the name of this transformation as a string.
- validate(node, options=None)[source]#
Ensure that it is valid to apply this transformation to the supplied node.
- Parameters:
node (
psyclone.psyir.nodes.Routine) – the node to validate.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- Raises:
TransformationError – if the node is not a Routine.
- class psyclone.psyir.transformations.HoistLocalArraysTrans[source]#
This transformation takes a Routine and promotes any local arrays to the Container scope:
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Assignment >>> from psyclone.psyir.transformations import HoistLocalArraysTrans >>> code = ("module test_mod\n" ... "contains\n" ... " subroutine test_sub(n)\n" ... " integer :: i,j,n\n" ... " real :: a(n,n)\n" ... " real :: value = 1.0\n" ... " do i=1,n\n" ... " do j=1,n\n" ... " a(i,j) = value\n" ... " end do\n" ... " end do\n" ... " end subroutine test_sub\n" ... "end module test_mod\n") >>> psyir = FortranReader().psyir_from_source(code) >>> hoist = HoistLocalArraysTrans() >>> hoist.apply(psyir.walk(Routine)[0]) >>> print(FortranWriter()(psyir).lower()) module test_mod implicit none real, allocatable, dimension(:,:), private :: a public public :: test_sub contains subroutine test_sub(n) integer :: n integer :: i integer :: j real :: value = 1.0 if (.not.allocated(a) .or. ubound(a, 1) /= n .or. ubound(a, 2) /= n) then if (allocated(a)) then deallocate(a) end if allocate(a(1 : n, 1 : n)) end if do i = 1, n, 1 do j = 1, n, 1 a(i,j) = value enddo enddo end subroutine test_sub end module test_mod
By default, the target routine will be rejected if it is found to contain an ACCRoutineDirective since this usually implies that the routine will be launched in parallel on the OpenACC device. This check can be disabled by setting ‘allow_accroutine’ to True in the options dictionary.
Inheritance

- apply(node, options=None)[source]#
Applies the transformation to the supplied Routine node, moving any local arrays up to Container scope and adding a suitable allocation when they are first accessed. If there are no local arrays or the supplied Routine is a program then this method does nothing.
- Parameters:
node (
psyclone.psyir.nodes.Routine) – target PSyIR node.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
options["allow_accroutine"] (bool) – permit the target routine to contain an ACCRoutineDirective. These are forbidden by default because their presence usually indicates that the routine will be run in parallel on the OpenACC device.
- validate(node, options=None)[source]#
Checks that the supplied node is a valid target for a hoist- local-arrays transformation. It must be a Routine that is within a Container (that is not a FileContainer).
- Parameters:
node (subclass of
psyclone.psyir.nodes.Routine) – target PSyIR node.options (Optional[Dict[str, Any]]) – any options for the transformation.
- Raises:
TransformationError – if the supplied node is not a Routine.
TransformationError – if the Routine is not within a Container (that is not a FileContainer).
TransformationError – if the routine contains an OpenACC routine directive and options[‘allow_accroutine’] is not True.
TransformationError – if any symbols corresponding to local arrays have a tag that already exists in the table of the parent Container.
- class psyclone.psyir.transformations.HoistLoopBoundExprTrans[source]#
This transformation moves complex bounds expressions out of the loop construct and places them in integer scalar assignments before the loop.
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Loop >>> from psyclone.psyir.transformations import HoistTrans >>> code = ("program test\n" ... " use mymod, only: mytype\n" ... " integer :: i,j,n\n" ... " real :: a(n)\n" ... " do i=mytype%start, UBOUND(a,1)\n" ... " a(i) = 1.0\n" ... " end do\n" ... "end program\n") >>> psyir = FortranReader().psyir_from_source(code) >>> hoist = HoistLoopBoundExprTrans() >>> hoist.apply(psyir.walk(Loop)[0]) >>> print(FortranWriter()(psyir)) program test use mymod, only : mytype integer :: i integer :: j integer :: n real, dimension(n) :: a integer :: loop_bound integer :: loop_bound_1 loop_bound_1 = UBOUND(a, 1) loop_bound = mytype%start do i = loop_bound, loop_bound_1, 1 a(i) = 1.0 enddo end program test
Inheritance

- apply(node, options=None)[source]#
Move complex bounds expressions out of the given loop construct and place them in integer scalar assignments before the loop.
- Parameters:
node (
psyclone.psyir.nodes.Loop) – target PSyIR loop.options (Dict[str, Any]) – a dictionary with options for transformations.
- validate(node, options=None)[source]#
Checks that the supplied node is a valid target for the transformation.
- Parameters:
node (
psyclone.psyir.nodes.Loop) – target PSyIR loop.options (Dict[str, Any]) – a dictionary with options for transformations.
- Raises:
TransformationError – if the supplied node does not have an ancestor Routine.
TransformationError – if the supplied node is directly inside a Directive Schedule.
- class psyclone.psyir.transformations.HoistTrans[source]#
This transformation takes an assignment and moves it outside of its parent loop if it is valid to do so. If as a result the loop body becomes empty, the loop will be removed altogether. For example:
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Assignment >>> from psyclone.psyir.transformations import HoistTrans >>> code = ("program test\n" ... " integer :: i,j,n\n" ... " real :: a(n,n)\n" ... " real value\n" ... " do i=1,n\n" ... " value = 1.0\n" ... " do j=1,n\n" ... " a(i,j) = value\n" ... " end do\n" ... " end do\n" ... "end program\n") >>> psyir = FortranReader().psyir_from_source(code) >>> hoist = HoistTrans() >>> hoist.apply(psyir.walk(Assignment)[0]) >>> print(FortranWriter()(psyir)) program test integer :: i integer :: j integer :: n real, dimension(n,n) :: a real :: value value = 1.0 do i = 1, n, 1 do j = 1, n, 1 a(i,j) = value enddo enddo end program test
Inheritance

- apply(node, options=None)[source]#
Applies the hoist transformation to the supplied assignment node within a loop, moving the assignment outside of the loop if it is valid to do so. Issue #1445 will also look to extend this transformation to other types of node.
- Parameters:
node (subclass of
psyclone.psyir.nodes.Assignment) – target PSyIR node.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- validate(node, options=None)[source]#
Checks that the supplied node is a valid target for a hoist transformation. At this stage only an assignment statement is allowed to be hoisted, see #1445. It should also be tested if there is a directive outside of the loop, see #1446
- Parameters:
node (subclass of
psyclone.psyir.nodes.Assignment) – target PSyIR node.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- Raises:
TransformationError – if the supplied node is not an assignment.
TransformationError – if the assignment is not within a loop.
TransformationError – if the assignment is not a direct child of the the loop.
- class psyclone.psyir.transformations.IncreaseRankLoopArraysTrans[source]#
This transformation takes a loop and a list of arrays accessed inside the loop, and increases those arrays with an additional dimension with the size of the interation space. Then it indexes all accesses with the loop variable, so that each iteration accesses a unique location. Effectively making the sub-array private for each iteration of the loop. It also indexes assignments outside the loop to iterate over the whole new rank.
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Loop >>> from psyclone.psyir.transformations import IncreaseRankLoopArraysTrans >>> code = (""" ... program test ... integer :: N=10, M=10 ... integer :: i, j ... real, dimension(N) :: ztmp ... ! (Array) Assignments to the target variable outside the loop are ... ! permited, but any other use will not pass the validation ... ztmp = 0 ... ztmp(0) = 1 ... do i = -5, M + 3 ... do j = 1, N ... ztmp(j) = 1 ... end do ... do j = 1, N ... ztmp(j) = ztmp(j) + 1 ... end do ... end do ... end program ... """) >>> psyir = FortranReader().psyir_from_source(code) >>> irla = IncreaseRankLoopArraysTrans() >>> irla.apply(psyir.walk(Loop)[0], arrays=['ztmp']) >>> print(FortranWriter()(psyir)) program test integer, save :: n = 10 integer, save :: m = 10 integer :: i integer :: j real, dimension(n,-5:m + 3) :: ztmp ztmp = 0 ztmp(0,:) = 1 do i = -5, m + 3, 1 do j = 1, n, 1 ztmp(j,i) = 1 enddo do j = 1, n, 1 ztmp(j,i) = ztmp(j,i) + 1 enddo enddo end program test
Inheritance

- validate(node, **kwargs)[source]#
Checks that the supplied node is a valid target.
- Parameters:
node (
Loop) – target Loop node.arrays (list[psyclone.psyir.symbols.symbol.Symbol | str] | None) – list of arrays that will have their rank increased.
- Raises:
TransformationError – if the node is not a Loop.
TransformationError – if the node is not inside a Routine.
TransformationError – if the target node does not have static loop bounds (which the dimension declaration needs).
TransformationError – if no array names are provided.
TransformationError – the given array name or the symbol is not local or not an array.
TransformationError – if any of the arrays are referenced inside a CodeBlock.
TransformationError – if any of the arrays are referenced outside the loop in anything other than (array) assignments to their variable.
- class psyclone.psyir.transformations.InlineTrans[source]#
This transformation takes a Call (which may have a return value) and replaces it with the body of the target routine. It is used as follows:
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Call, Routine >>> from psyclone.psyir.transformations import InlineTrans >>> code = """ ... module test_mod ... contains ... subroutine run_it() ... integer :: i ... real :: a(10) ... do i=1,10 ... a(i) = 1.0 ... call sub(a(i)) ... end do ... end subroutine run_it ... subroutine sub(x) ... real, intent(inout) :: x ... x = 2.0*x ... end subroutine sub ... end module test_mod""" >>> psyir = FortranReader().psyir_from_source(code) >>> call = psyir.walk(Call)[0] >>> inline_trans = InlineTrans() >>> inline_trans.apply(call) >>> # Uncomment the following line to see a text view of the schedule >>> # print(psyir.walk(Routine)[0].view()) >>> print(FortranWriter()(psyir.walk(Routine)[0])) subroutine run_it() integer :: i real, dimension(10) :: a do i = 1, 10, 1 a(i) = 1.0 a(i) = 2.0 * a(i) enddo end subroutine run_it
The target of the call must already be present in the same Container (module) as the call site. This may be achieved using KernelModuleInlineTrans.
Warning
Routines/calls with any of the following characteristics are not supported and will result in a TransformationError:
the routine contains an early Return statement;
the routine contains a variable with UnknownInterface;
the routine contains a variable with StaticInterface;
the routine contains an UnsupportedType variable with ArgumentInterface;
the shape of any array arguments as declared inside the routine does not match the shape of the arrays being passed as arguments;
the routine accesses an un-resolved symbol;
the routine accesses a symbol declared in the Container to which it belongs.
Some of these restrictions will be lifted by #924.
Inheritance

- apply(node, routine=None, use_first_callee_and_no_arg_check=False, permit_codeblocks=False, permit_unsupported_type_args=False, **kwargs)[source]#
Takes the body of the routine that is the target of the supplied call and replaces the call with it.
- Parameters:
node (
Call) – the Call node to inline.routine (
Optional[Routine]) – Optional Routine to be inlined. (By default, PSyclone will search for a target routine with a matching signature).use_first_callee_and_no_arg_check (
bool) – if True, simply use the first potential callee routine. No argument type-checking is performed.permit_codeblocks (
bool) – If False (the default), raise an Exception if the target routine contains a CodeBlock.permit_unsupported_type_args (
bool) – If True then the target routine is permitted to have arguments of UnsupportedType.
- Raises:
InternalError – if the merge of the symbol tables fails. In theory this should never happen because validate() should catch such a situation.
- validate(node, **kwargs)[source]#
Checks that the supplied node is a valid target for inlining.
- Parameters:
node (
Call) – the target PSyIR Call to inline.routine (psyclone.psyir.nodes.routine.Routine | None) – Optional Routine to be inlined. (By default, PSyclone will search for a target routine with a matching signature).
use_first_callee_and_no_arg_check (bool) – if True, simply use the first potential callee routine. No argument type-checking is performed.
permit_codeblocks (bool) – If False (the default), raise an Exception if the target routine contains a CodeBlock.
permit_unsupported_type_args (bool) – If True then the target routine is permitted to have arguments of UnsupportedType.
- Raises:
TransformationError – if the supplied node is not a Call, is an IntrinsicCall or a call to a PSyclone-generated or ELEMENTAL routine.
TransformationError – if the routine body contains a Return that is not the first or last statement.
TransformationError – if the routine body contains a CodeBlock and permit_codeblocks is not True.
TransformationError – if the routine contains an ALLOCATE (as we don’t support adding the required DEALLOCATE).
TransformationError – if the called routine has a named argument.
TransformationError – if the call-site is not within a Routine.
TransformationError – if any of the variables declared within the called routine are of UnknownInterface.
TransformationError – if any of the variables declared within the called routine have a StaticInterface.
TransformationError – if any of the subroutine arguments is of UnsupportedType.
TransformationError – if a symbol of a given name is imported from different containers at the call site and within the routine.
TransformationError – if the routine accesses an un-resolved symbol.
TransformationError – if the number of arguments in the call does not match the number of formal arguments of the routine.
TransformationError – if a symbol declared in the parent container is accessed in the target routine.
TransformationError – if the shape of an array formal argument does not match that of the corresponding actual argument.
TransformationError – if one of the declaratoins in the routine depends on an argument that is written to prior to the call.
InternalError – if an unhandled Node type is returned by Reference.previous_accesses().
- Return type:
None
- class psyclone.psyir.transformations.Abs2CodeTrans[source]#
Provides a transformation from a PSyIR ABS Operator node to equivalent code in a PSyIR tree. Validity checks are also performed.
The transformation replaces
R = ABS(X)
with the following logic:
IF X < 0.0: R = X*-1.0 ELSE: R = X
Inheritance

- apply(node, options=None, **kwargs)[source]#
Apply the ABS intrinsic conversion transformation to the specified node. This node must be an ABS UnaryOperation. The ABS UnaryOperation is converted to equivalent inline code. This is implemented as a PSyIR transform from:
R = ... ABS(X) ...
to:
tmp_abs = X if tmp_abs < 0.0: res_abs = tmp_abs*-1.0 else: res_abs = tmp_abs R = ... res_abs ...
where
Xcould be an arbitrarily complex PSyIR expression and...could be arbitrary PSyIR code.This transformation requires the operation node to be a descendant of an assignment and will raise an exception if this is not the case.
- Parameters:
node (
psyclone.psyir.nodes.UnaryOperation) – an ABS UnaryOperation node.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- validate(node, options=None)[source]#
Check that it is safe to apply the transformation to the supplied node.
- Parameters:
node (
psyclone.psyir.nodes.IntrinsicCall) – the SIGN call to transform.options (dict[str, Any]) – any of options for the transformation.
- class psyclone.psyir.transformations.DotProduct2CodeTrans[source]#
Provides a transformation from a PSyIR DOT_PRODUCT Operator node to equivalent code in a PSyIR tree. Validity checks are also performed.
If
Ris a scalar andA, andBhave dimensionN, The transformation replaces:R = ... DOT_PRODUCT(A,B) ...
with the following code:
TMP = 0.0 do I=1,N TMP = TMP + A(i)*B(i) R = ... TMP ...
For example:
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import IntrinsicCall >>> from psyclone.psyir.transformations import DotProduct2CodeTrans >>> code = ("subroutine dot_product_test(v1,v2)\n" ... "real,intent(in) :: v1(10), v2(10)\n" ... "real :: result\n" ... "result = dot_product(v1,v2)\n" ... "end subroutine\n") >>> psyir = FortranReader().psyir_from_source(code) >>> trans = DotProduct2CodeTrans() >>> trans.apply(psyir.walk(IntrinsicCall)[0]) >>> print(FortranWriter()(psyir)) subroutine dot_product_test(v1, v2) real, dimension(10), intent(in) :: v1 real, dimension(10), intent(in) :: v2 real :: result integer :: i real :: res_dot_product res_dot_product = 0.0 do i = 1, 10, 1 res_dot_product = res_dot_product + v1(i) * v2(i) enddo result = res_dot_product end subroutine dot_product_test
Inheritance

- apply(node, options=None, **kwargs)[source]#
Apply the DOT_PRODUCT intrinsic conversion transformation to the specified node. This node must be a DOT_PRODUCT BinaryOperation. If the transformation is successful then an assignment which includes a DOT_PRODUCT BinaryOperation node is converted to equivalent inline code.
- Parameters:
node (
psyclone.psyir.nodes.BinaryOperation) – a DOT_PRODUCT Binary-Operation node.options (dict of str:str or None) – a dictionary with options for transformations.
- validate(node, options=None, **kwargs)[source]#
Perform checks to ensure that it is valid to apply the DotProduct2CodeTran transformation to the supplied node.
Note, this validation does not check for invalid argument combinations to dot_product (e.g. different precision or different datatypes) as that should have already been picked up when creating the PSyIR.
- Parameters:
node (
psyclone.psyir.nodes.BinaryOperation) – the node that is being checked.options (dict of str:str or None) – a dictionary with options for transformations.
- Raises:
TransformationError – if one of the arguments is not a Reference node.
TransformationError – if an argument does not use array slice notation and is not a 1d array.
TransformationError – if an argument uses array slice notation but the array slice is not for the first dimension of the array.
TransformationError – if an argument uses array slice notation but it is not for the full range of the dimension.
- class psyclone.psyir.transformations.Matmul2CodeTrans[source]#
Provides a transformation from a PSyIR MATMUL Operator node to equivalent code in a PSyIR tree. Validity checks are also performed.
For a matrix-vector multiplication, if the dimensions of
R,A, andBareR(N),A(N,M),B(M), the transformation replaces:R=MATMUL(A,B)
with the following code:
do i=1,N R(i) = 0.0 do j=1,M R(i) = R(i) + A(i,j) * B(j)
For a matrix-matrix multiplication, if the dimensions of
R,A, andBareR(P,M),A(P,N),B(N,M), the MATMUL is replaced with the following code:do j=1,M do i=1,P R(i,j) = 0.0 do ii=1,N R(i,j) = R(i,j) + A(i,ii) * B(ii,j)
Note that this transformation does not support the case where
Ais a rank-1 array.Inheritance

- apply(node, options=None, **kwargs)[source]#
Apply the MATMUL intrinsic conversion transformation to the specified node. This node must be a MATMUL IntrinsicCall. The first argument must currently have two dimensions while the second must have either one or two dimensions. Each argument is permitted to have additional dimensions (i.e. more than 2) but in each case it is only the first one or two which may be ranges. Further, the ranges must currently be for the full index space for that dimension (i.e. array subsections are not supported). If the transformation is successful then an assignment which includes a MATMUL IntrinsicCall node is converted to equivalent inline code.
- Parameters:
node (
psyclone.psyir.nodes.IntrinsicCall) – a MATMUL IntrinsicCall node.options (Optional[Dict[str, Any]]) – options for the transformation.
- validate(node, options=None, **kwargs)[source]#
Perform checks to ensure that it is valid to apply the Matmul2CodeTran transformation to the supplied node.
- Parameters:
node (
psyclone.psyir.nodes.IntrinsicCall) – the node that is being checked.options (Optional[Dict[str, Any]]) – options for the transformation.
- Raises:
TransformationError – if the node argument is not the expected type.
TransformationError – if the parent of the MATMUL operation is not an assignment.
TransformationError – if the matmul arguments are not in the required form.
TransformationError – if sub-sections of an array are present in the arguments.
- class psyclone.psyir.transformations.Max2CodeTrans[source]#
Provides a transformation from a PSyIR MAX Intrinsic node to equivalent code in a PSyIR tree. Validity checks are also performed (by a parent class).
The transformation replaces
R = MAX(A, B, C ...)
with the following logic:
R = A if B > R: R = B if C > R: R = C ...
Inheritance

- apply(node, options=None, **kwargs)[source]#
Applies the Max2CodeTrans to the provided node.
- Parameters:
node (
psyclone.psyir.nodes.IntrinsicCall) – a MAX intrinsic.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- class psyclone.psyir.transformations.Min2CodeTrans[source]#
Provides a transformation from a PSyIR MIN Intrinsic node to equivalent code in a PSyIR tree. Validity checks are also performed (by a parent class).
The transformation replaces
R = MIN(A, B, C ...)
with the following logic:
R = A if B < R: R = B if C < R: R = C ...
Inheritance

- apply(node, options=None, **kwargs)[source]#
Applies the Min2CodeTrans to the provided node.
- Parameters:
node (
psyclone.psyir.nodes.IntrinsicCall) – a MIN intrinsic.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- class psyclone.psyir.transformations.Sign2CodeTrans[source]#
Provides a transformation from a PSyIR SIGN intrinsic node to equivalent code in a PSyIR tree. Validity checks are also performed.
The transformation replaces
R = SIGN(A, B)
with the following logic:
R = ABS(A) if B < 0.0: R = R*-1.0
i.e. the value of
Awith the sign ofBInheritance

- apply(node, options=None, **kwargs)[source]#
Apply the SIGN intrinsic conversion transformation to the specified node. This node must be a SIGN IntrinsicCall. The SIGN IntrinsicCall is converted to equivalent inline code. This is implemented as a PSyIR transform from:
R = ... SIGN(A, B) ...
to:
tmp_abs = A if tmp_abs < 0.0: res_abs = tmp_abs*-1.0 else: res_abs = tmp_abs res_sign = res_abs tmp_sign = B if tmp_sign < 0.0: res_sign = res_sign*-1.0 R = ... res_sign ...
where
AandBcould be arbitrarily complex PSyIR expressions,...could be arbitrary PSyIR code and whereABShas been replaced with inline code by the NemoAbsTrans transformation.This transformation requires the IntrinsicCall node to be a child of an assignment and will raise an exception if this is not the case.
- Parameters:
node (
psyclone.psyir.nodes.IntrinsicCall) – a SIGN IntrinsicCall node.symbol_table (
psyclone.psyir.symbols.SymbolTable) – the symbol table.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- validate(node, options=None, **kwargs)[source]#
Check that it is safe to apply the transformation to the supplied node.
- Parameters:
node (
psyclone.psyir.nodes.IntrinsicCall) – the SIGN call to transform.options (dict[str, Any]) – any options for the transformation.
symbol_table (
psyclone.psyir.symbols.SymbolTable) – the symbol table.
- class psyclone.psyir.transformations.Sum2LoopTrans[source]#
Provides a transformation from a PSyIR SUM IntrinsicCall node to an equivalent PSyIR loop structure that is suitable for running in parallel on CPUs and GPUs. Validity checks are also performed.
If SUM contains a single positional argument which is an array, all elements of that array are summed and the result returned in the scalar R.
R = SUM(ARRAY)
For example, if the array is two dimensional, the equivalent code for real data is:
R = 0.0 DO J=LBOUND(ARRAY,2),UBOUND(ARRAY,2) DO I=LBOUND(ARRAY,1),UBOUND(ARRAY,1) R = R + ARRAY(I,J)
If the mask argument is provided then the mask is used to determine whether the sum is applied:
R = SUM(ARRAY, mask=MOD(ARRAY, 2.0)==1)
If the array is two dimensional, the equivalent code for real data is:
R = 0.0 DO J=LBOUND(ARRAY,2),UBOUND(ARRAY,2) DO I=LBOUND(ARRAY,1),UBOUND(ARRAY,1) IF (MOD(ARRAY(I,J), 2.0)==1) THEN R = R + ARRAY(I,J)
The dimension argument is currently not supported and will result in a TransformationError exception being raised.
R = SUM(ARRAY, dimension=2)
The array passed to MAXVAL may use any combination of array syntax, array notation, array sections and scalar bounds:
R = SUM(ARRAY) ! array syntax R = SUM(ARRAY(:,:)) ! array notation R = SUM(ARRAY(1:10,lo:hi)) ! array sections R = SUM(ARRAY(1:10,:)) ! mix of array section and array notation R = SUM(ARRAY(1:10,2)) ! mix of array section and scalar bound
For example:
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.transformations import Sum2LoopTrans >>> code = ("subroutine sum_test(array,n,m)\n" ... " integer :: n, m\n" ... " real :: array(10,10)\n" ... " real :: result\n" ... " result = sum(array)\n" ... "end subroutine\n") >>> psyir = FortranReader().psyir_from_source(code) >>> sum_node = psyir.children[0].children[0].children[1] >>> Sum2LoopTrans().apply(sum_node) >>> print(FortranWriter()(psyir)) subroutine sum_test(array, n, m) integer :: n integer :: m real, dimension(10,10) :: array real :: result integer :: idx integer :: idx_1 result = 0.0 do idx = 1, 10, 1 do idx_1 = 1, 10, 1 result = result + array(idx_1,idx) enddo enddo end subroutine sum_test
Inheritance

- apply(node, options=None, **kwargs)[source]#
Apply the Sum2LoopTrans to the input node.
- Parameters:
node (
psyclone.psyir.nodes.IntrinsicCall) – the SUM intrinsic to transform.options (Optional[Dict[str, Any]]) – options for the transformation.
- class psyclone.psyir.transformations.LoopFuseTrans[source]#
Provides a generic loop-fuse transformation to two Nodes in the PSyIR of a Schedule after performing validity checks for the supplied Nodes. Examples are given in the descriptions of any children classes.
If loops have different named loop variables, when possible the loop variable of the second loop will be renamed to be the same as the first loop. This has the side effect that the second loop’s variable will no longer have its value modified, with the expectation that that value isn’t used anymore.
Note that the validation of this transformation still has several shortcomings, especially for domain API loops. Use at your own risk.
Inheritance

- apply(node1, node2, options=None)[source]#
Fuses two loops represented by psyclone.psyir.nodes.Node objects after performing validity checks.
If the two loops don’t have the same loop variable, the second loop’s variable (and references to it inside the loop) will be changed to be references to the first loop’s variable before merging. This has the side effect that the second loop’s variable will no longer have its value modified, with the expectation that that value isn’t used after.
- Parameters:
node1 (
psyclone.psyir.nodes.Node) – the first Node that is being checked.node2 (
psyclone.psyir.nodes.Node) – the second Node that is being checked.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- validate(node1, node2, options=None)[source]#
Performs various checks to ensure that it is valid to apply the LoopFuseTrans transformation to the supplied Nodes.
- Parameters:
node1 (
psyclone.psyir.nodes.Node) – the first Node that is being checked.node2 (
psyclone.psyir.nodes.Node) – the second Node that is being checked.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
options["force"] (bool) – whether to force fusion of the target loop (i.e. ignore any dependence analysis). This only skips a limited number of the checks, and does not fully force merging.
- Raises:
TransformationError – if one or both of the Nodes is/are not a
psyclone.psyir.nodes.Loop.TransformationError – if one or both Nodes are not fully-formed.
TransformationError – if the Nodes do not have the same parent.
TransformationError – if the Nodes are not next to each other in the tree.
TransformationError – if the two Loops do not have the same iteration space.
TransformationError – if the two Loops don’t share the same iteration variable, but make use of the other loop’s variable in their body.
TransformationError – if there are dependencies between the loops that prevent the loop fusion.
- class psyclone.psyir.transformations.LoopSwapTrans[source]#
Provides a loop-swap transformation, e.g.:
DO j=1, m DO i=1, n
becomes:
DO i=1, n DO j=1, m
This transform is used as follows:
>>> from psyclone.parse.algorithm import parse >>> from psyclone.psyGen import PSyFactory >>> ast, invokeInfo = parse("shallow_alg.f90") >>> psy = PSyFactory("gocean").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 LoopSwapTrans >>> swap = LoopSwapTrans() >>> swap.apply(schedule.children[0]) >>> # Uncomment the following line to see a text view of the schedule >>> # print(schedule.view())
Inheritance

- apply(node, options=None)[source]#
The argument
outermust be a loop which has exactly one inner loop. This transform then swaps the outer and inner loop.- Parameters:
outer (
psyclone.psyir.nodes.Loop) – the node representing the outer loop.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- Raises:
TransformationError – if the supplied node does not allow a loop swap to be done.
- validate(node, options=None)[source]#
Checks if the given node contains a valid Fortran structure to allow swapping loops. This means the node must represent a loop, and it must have exactly one child that is also a loop.
- Parameters:
node_outer (py:class:psyclone.psyir.nodes.Loop) – a Loop node from an AST.
options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- Raises:
TransformationError – if the supplied node does not allow a loop swap to be done.
TransformationError – if either the inner or outer loop has a symbol table.
- class psyclone.psyir.transformations.LoopTilingTrans[source]#
Apply a loop tiling transformation to a loop. For example:
>>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Loop >>> from psyclone.psyir.transformations import LoopTilingTrans >>> psyir = FortranReader().psyir_from_source(""" ... subroutine sub() ... integer :: i, j, tmp(100) ... do i=1, 100 ... do j=1, 100 ... tmp(i, j) = 2 * tmp(i, j) ... enddo ... enddo ... end subroutine sub""") >>> loop = psyir.walk(Loop)[0] >>> LoopTilingTrans().apply(loop, tiledims=[32, 32])
will generate:
subroutine sub() integer :: i integer :: j integer, dimension(100) :: tmp integer :: j_out_var integer :: i_out_var do i_out_var = 1, 100, 32 do j_out_var = 1, 100, 32 do i = i_out_var, MIN(i_out_var + (32 - 1), 100), 1 do j = j_out_var, MIN(j_out_var + (32 - 1), 100), 1 tmp(i, j) = 2 * tmp(i, j) enddo enddo enddo enddo end subroutine sub
Inheritance

- apply(outer_loop, tiledims, **kwargs)[source]#
Converts the given Loop construct into a tiled version of the nested loops.
- Parameters:
outer_loop (
Loop) – the loop to transform.tiledims (
list[int]) – the dimensions of the tile.
- Raises:
TransformationError – if any of the tile dimensions are not positive integers
TransformationError – if the transformation cannot be applied
- validate(outer_loop, tiledims, **kwargs)[source]#
Validates that the given Loop node can have a LoopTilingTrans applied.
- Parameters:
outer_loop (
Loop) – the loop to validate.tiledims (
list[int]) – the dimensions of the tile.
- Raises:
TransformationError – if any of the tile dimensions are not positive integers
TransformationError – if the transformation cannot be applied
- class psyclone.psyir.transformations.LoopTiling2DTrans[source]#
Apply a 2D loop tiling transformation to a loop. This is a special case of LoopTilingTrans for 2D square tiles. For example:
>>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Loop >>> from psyclone.psyir.transformations import LoopTiling2DTrans >>> psyir = FortranReader().psyir_from_source(""" ... subroutine sub() ... integer :: i, j, tmp(100) ... do i=1, 100 ... do j=1, 100 ... tmp(i, j) = 2 * tmp(i, j) ... enddo ... enddo ... end subroutine sub""") >>> loop = psyir.walk(Loop)[0] >>> LoopTiling2DTrans().apply(loop)
will generate:
subroutine sub() integer :: i integer :: j integer, dimension(100) :: tmp integer :: j_out_var integer :: i_out_var do i_out_var = 1, 100, 32 do j_out_var = 1, 100, 32 do i = i_out_var, MIN(i_out_var + (32 - 1), 100), 1 do j = j_out_var, MIN(j_out_var + (32 - 1), 100), 1 tmp(i, j) = 2 * tmp(i, j) enddo enddo enddo enddo end subroutine sub
Inheritance

- apply(node, options=None)[source]#
Converts the given 2D Loop construct into a tiled version of the nested loops.
- Parameters:
node (
psyclone.psyir.nodes.Loop) – the loop to transform.options (Optional[Dict[str, Any]]) – a dict with options for transformations.
options["tilesize"] (int) – The size of the resulting tile, currently square tiles are always used. If not specified, the value 32 is used.
- validate(node, options=None)[source]#
Validates that the given Loop node can have a LoopTiling2DTrans applied.
- Parameters:
node (
psyclone.psyir.nodes.Loop) – the loop to validate.options (Optional[Dict[str, Any]]) – a dict with options for transformation.
options["tilesize"] (int) – The size of the resulting tile, currently square tiles are always used. If not specified, the value 32 is used.
- Raises:
TransformationError – if an unsupported option has been provided.
TransformationError – if the provided tilesize is not a integer.
- class psyclone.psyir.transformations.LoopTrans[source]#
This abstract class is a base class for all transformations that act on a Loop node. It gives access to a validate function that makes sure that the supplied node is a Loop. We also check that all nodes to be enclosed are valid for this transformation - this requires that the sub-class populate the excluded_node_types tuple.
Inheritance

- apply(node, options=None, node_type_check=True, verbose=False, **kwargs)[source]#
Applies the transformation to the provided node.
This function only calls the superclass method, but is required for option specification.
- Parameters:
node (subclass of
psyclone.psyir.nodes.Node) – target PSyIR node.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
node_type_check (
bool) – If the type of nodes enclosed in the loop should be tested to avoid including unsupported nodes in the transformation.verbose (
bool) – whether to log the reason the validation failed, at the moment with a comment in the provided PSyIR node.
- property name#
- Returns:
the name of this class.
- Return type:
str
- validate(node, options=None, **kwargs)[source]#
Checks that the supplied node is a valid target for a loop transformation.
- Parameters:
node (subclass of
psyclone.psyir.nodes.Node) – target PSyIR node.- Raises:
TransformationError – if the supplied node is not a (fully- formed) Loop.
TransformationError – if any of the nodes within the loop are of an unsupported type.
TransformationError – if the loop is of ‘null’ type.
TransformationError – if the supplied options are not a dictionary.
- class psyclone.psyir.transformations.Maxval2LoopTrans[source]#
Provides a transformation from a PSyIR MAXVAL IntrinsicCall node to an equivalent PSyIR loop structure that is suitable for running in parallel on CPUs and GPUs. Validity checks are also performed.
If MAXVAL contains a single positional argument which is an array, the maximum value of all of the elements in the array is returned in the the scalar R.
R = MAXVAL(ARRAY)
For example, if the array is two dimensional, the equivalent code for real data is:
R = -HUGE(R) DO J=LBOUND(ARRAY,2),UBOUND(ARRAY,2) DO I=LBOUND(ARRAY,1),UBOUND(ARRAY,1) R = MAX(R, ARRAY(I,J))
If the mask argument is provided then the mask is used to determine whether the maxval is applied:
R = MAXVAL(ARRAY, mask=MOD(ARRAY, 2.0)==1)
If the array is two dimensional, the equivalent code for real data is:
R = -HUGE(R) DO J=LBOUND(ARRAY,2),UBOUND(ARRAY,2) DO I=LBOUND(ARRAY,1),UBOUND(ARRAY,1) IF (MOD(ARRAY(I,J), 2.0)==1) THEN R = MAX(R, ARRAY(I,J))
The dimension argument is currently not supported and will result in a TransformationError exception being raised.
R = MAXVAL(ARRAY, dimension=2)
The array passed to MAXVAL may use any combination of array syntax, array notation, array sections and scalar bounds:
R = MAXVAL(ARRAY) ! array syntax R = MAXVAL(ARRAY(:,:)) ! array notation R = MAXVAL(ARRAY(1:10,lo:hi)) ! array sections R = MAXVAL(ARRAY(1:10,:)) ! mix of array section and array notation R = MAXVAL(ARRAY(1:10,2)) ! mix of array section and scalar bound
An example use of this transformation is given below:
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.transformations import Maxval2LoopTrans >>> code = ("subroutine maxval_test(array)\n" ... " real :: array(10,10)\n" ... " real :: result\n" ... " result = maxval(array)\n" ... "end subroutine\n") >>> psyir = FortranReader().psyir_from_source(code) >>> sum_node = psyir.children[0].children[0].children[1] >>> Maxval2LoopTrans().apply(sum_node) >>> print(FortranWriter()(psyir)) subroutine maxval_test(array) real, dimension(10,10) :: array real :: result integer :: idx integer :: idx_1 result = -HUGE(result) do idx = 1, 10, 1 do idx_1 = 1, 10, 1 result = MAX(result, array(idx_1,idx)) enddo enddo end subroutine maxval_test
Inheritance

- apply(node, options=None, **kwargs)[source]#
Apply the Maxval2LoopTrans to the input node.
- Parameters:
node (
psyclone.psyir.nodes.IntrinsicCall) – the MAXVAL intrinsic to transform.options (Optional[Dict[str, Any]]) – options for the transformation.
- class psyclone.psyir.transformations.Minval2LoopTrans[source]#
Provides a transformation from a PSyIR MINVAL IntrinsicCall node to an equivalent PSyIR loop structure that is suitable for running in parallel on CPUs and GPUs. Validity checks are also performed.
If MINVAL contains a single positional argument which is an array, the minimum value of all of the elements in the array is returned in the the scalar R.
R = MINVAL(ARRAY)
For example, if the array is two dimensional, the equivalent code for real data is:
R = HUGE(R) DO J=LBOUND(ARRAY,2),UBOUND(ARRAY,2) DO I=LBOUND(ARRAY,1),UBOUND(ARRAY,1) R = MIN(R, ARRAY(I,J))
If the mask argument is provided then the mask is used to determine whether the minval is applied:
R = MINVAL(ARRAY, mask=MOD(ARRAY, 2.0)==1)
If the array is two dimensional, the equivalent code for real data is:
R = HUGE(R) DO J=LBOUND(ARRAY,2),UBOUND(ARRAY,2) DO I=LBOUND(ARRAY,1),UBOUND(ARRAY,1) IF (MOD(ARRAY(I,J), 2.0)==1) THEN R = MIN(R, ARRAY(I,J))
The dimension argument is currently not supported and will result in a TransformationError exception being raised.
R = MINVAL(ARRAY, dimension=2)
The array passed to MINVAL may use any combination of array syntax, array notation, array sections and scalar bounds:
R = MINVAL(ARRAY) ! array syntax R = MINVAL(ARRAY(:,:)) ! array notation R = MINVAL(ARRAY(1:10,lo:hi)) ! array sections R = MINVAL(ARRAY(1:10,:)) ! mix of array section and array notation R = MINVAL(ARRAY(1:10,2)) ! mix of array section and scalar bound
For example:
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.transformations import Minval2LoopTrans >>> code = ("subroutine minval_test(array)\n" ... " real :: array(10,10)\n" ... " real :: result\n" ... " result = minval(array)\n" ... "end subroutine\n") >>> psyir = FortranReader().psyir_from_source(code) >>> sum_node = psyir.children[0].children[0].children[1] >>> Minval2LoopTrans().apply(sum_node) >>> print(FortranWriter()(psyir)) subroutine minval_test(array) real, dimension(10,10) :: array real :: result integer :: idx integer :: idx_1 result = HUGE(result) do idx = 1, 10, 1 do idx_1 = 1, 10, 1 result = MIN(result, array(idx_1,idx)) enddo enddo end subroutine minval_test
Inheritance

- apply(node, options=None, **kwargs)[source]#
Apply the Minval2LoopTrans to the input node.
- Parameters:
node (
psyclone.psyir.nodes.IntrinsicCall) – the MINVAL intrinsic to transform.options (Optional[Dict[str, Any]]) – options for the transformation.
- class psyclone.psyir.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 “do”. The available options are “do” for “omp do”; “paralleldo” for “omp parallel do”; “teamsdistributeparalleldo” for “omp teams distribute parallel do”; “teamsloop” for “omp teams loop”; and “loop” for “omp loop”.
For example:
>>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.nodes import Loop >>> from psyclone.psyir.transformations import OMPLoopTrans >>> from psyclone.transformations import 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

- apply(node, options=None, reprod=None, enable_reductions=False, **kwargs)[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 transformationreprod (
bool) – indicating whether reproducible reductions should be used. By default the value from the config file will be used.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations and validation.
enable_reductions (
bool) – whether to attempt to infer reduction clauses or not.collapse (int | bool) – if it’s a bool and is False (default), it won’t collapse. If it’s a bool and is True, it will collapse as much as possible. If it’s an integer, it will attempt to collapse until the specified number of loops (if they exist and are safe to collapse them). The options ‘ignore_dependencies_for’ and ‘force’ also affect the collapse applicability analysis.
force (bool) – whether to force parallelisation of the target loop (i.e. ignore any dependence analysis).
ignore_dependencies_for (None | List[str]) – whether to ignore some symbol names from the dependence analysis checks.
sequential (bool) – whether this is a sequential loop.
verbose (bool) – whether to report the reasons the validate and collapse steps have failed.
nowait (bool) – whether to add a nowait clause and a corresponding barrier (or equivalent) to enable asynchronous execution.
privatise_arrays (bool) – whether to make the write after write dependency symbols declared as private.
reduction_ops (List[psyclone.psyir.nodes.operation.Operator | psyclone.psyir.nodes.intrinsic_call.IntrinsicCall.Intrinsic]) – if non-empty, attempt parallelisation of loops by inferring reduction clauses involving any of the reduction operators in the list.
node_type_check (bool) – If the type of nodes enclosed in the loop should be tested to avoid including unsupported nodes in the transformation.
- 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.psyir.transformations.OMPMinimiseSyncTrans[source]#
Attempts to remove OMPTaskwaitDirective or OMPBarrierDirective nodes from a supplied region as long as any dependencies continue to be satisfied. The goal is to reduce the number of synchronisation points to enable better asynchronous execution of loops or kernels.
For example:
>>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.nodes import Loop >>> from psyclone.psyir.transformations import OMPLoopTrans >>> from psyclone.psyir.transformations import OMPMinimiseSyncTrans >>> from psyclone.transformations import OMPParallelTrans >>> >>> psyir = FortranReader().psyir_from_source(""" ... subroutine test ... integer, dimension(100) :: a,b ... integer :: i ... ... do i = 1, 100 ... a(i) = i ... end do ... ... do i = 1, 100 ... b(i) = i ... end do ... ... do i = 1, 100 ... b(i) = b(i) + 1 ... end do ... ... do i = 1, 100 ... a(i) = a(i) + 1 ... end do ... end subroutine ... """) >>> omplooptrans1 = OMPLoopTrans() >>> for loop in psyir.walk(Loop): ... omplooptrans1.apply(loop, nowait=True) >>> partrans = OMPParallelTrans() >>> partrans.apply(psyir.children[0].children[:]) >>> rbartrans = OMPMinimiseSyncTrans() >>> rbartrans.apply(psyir.children[0]) >>> print(FortranWriter()(psyir)) subroutine test() integer, dimension(100) :: a integer, dimension(100) :: b integer :: i !$omp parallel default(shared), private(i) !$omp do schedule(auto) do i = 1, 100, 1 a(i) = i enddo !$omp end do nowait !$omp do schedule(auto) do i = 1, 100, 1 b(i) = i enddo !$omp end do nowait !$omp barrier !$omp do schedule(auto) do i = 1, 100, 1 b(i) = b(i) + 1 enddo !$omp end do nowait !$omp do schedule(auto) do i = 1, 100, 1 a(i) = a(i) + 1 enddo !$omp end do nowait !$omp barrier !$omp end parallel end subroutine test
Inheritance

- class psyclone.psyir.transformations.OMPTargetTrans[source]#
Adds an OpenMP target directive to a region of code.
For example:
>>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.nodes import Loop >>> from psyclone.psyir.transformations import OMPTargetTrans >>> >>> tree = 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 ... """) >>> OMPTargetTrans().apply(tree.walk(Loop)[0]) >>> print(FortranWriter()(tree)) subroutine my_subroutine() integer, dimension(10,10) :: a integer :: i integer :: j !$omp target do i = 1, 10, 1 do j = 1, 10, 1 a(i,j) = 0 enddo enddo !$omp end target end subroutine my_subroutine
Inheritance

- apply(node, options=None)[source]#
Insert an OMPTargetDirective before the provided node or list of nodes.
- Parameters:
node (List[
psyclone.psyir.nodes.Node]) – the PSyIR node or nodes to enclose in the OpenMP target region.options (Optional[Dict[str,Any]]) – a dictionary with options for transformations.
options["nowait"] (bool) – whether to add a nowait clause and a corresponding barrier to enable asynchronous execution.
options["device_string"] (str) – provide a compiler-platform identifier.
- validate(node, options=None)[source]#
Check that we can safely enclose the supplied node or list of nodes within an OpenMPTargetDirective.
- Parameters:
node (List[
psyclone.psyir.nodes.Node]) – the PSyIR node or nodes to enclose in the OpenMP target region.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
options["device_string"] (str) – provide a compiler-platform identifier.
options["allow_strings"] (str) – permit OMP target regions enclosing string operations.
options["verbose"] (str) – insert preceding comments with the reason that made this validation fail.
- Raises:
TransformationError – if it contains calls to routines that are not available in the accelerator device.
TransformationError – if its a function and the target region attempts to enclose the assingment setting the return value.
TransformationError – if the target region attempts to enclose string operations and the ‘allow_strings’ option is not set.
- class psyclone.psyir.transformations.OMPTaskTrans[source]#
Apply an OpenMP Task Transformation to a Loop. The Loop must be within an OpenMP Serial region (Single or Master) at codegen time. Once lowering begins, no more modifications to the tree should occur as the task directives do not recompute dependencies after lowering. In the future it may be possible to do this through an _update_node implementation.
Inheritance

- apply(node, options=None, **kwargs)[source]#
Apply the OMPTaskTrans to the specified node in a Schedule.
Can only be applied to a Loop.
The specified node is wrapped by directives during code generation like so:
!$OMP TASK ... !$OMP END TASK
At code-generation time, this node must be within (i.e. a child of) an OpenMP Serial region (OpenMP Single or OpenMP Master)
Any kernels or Calls will be inlined into the region before the task transformation is applied.
- Parameters:
node (
psyclone.psyir.nodes.Loop) – the supplied node to which we will apply the OMPTaskTrans transformationoptions (dictionary of string:values or None) – a dictionary with options for transformations and validation.
- property name#
- Returns:
the name of this transformation.
- Return type:
str
- validate(node, options=None, **kwargs)[source]#
Validity checks for input arguments.
Note that currently this implementation calls the apply methods of KernelModuleInlineTrans and FoldConditionalReturnExpressionsTrans. Although a copy of the provided node is used for this (so as to avoid any modifications to it in case the validation fails), this pattern indicates we should possibly re-write this transformation to require that those transformations have been applied first.
- Parameters:
node (
psyclone.psyir.nodes.Loop) – the Loop node to validate.options (dict of string:values or None) – a dictionary with options for transformations.
- class psyclone.psyir.transformations.OMPTaskwaitTrans[source]#
Adds zero or more OpenMP Taskwait directives to an OMP parallel region. This transformation will add directives to satisfy dependencies between Taskloop directives without an associated taskgroup (i.e. no nogroup clause). It also tries to minimise the number added to maximise available parallelism.
For example:
>>> from pysclone.parse.algorithm import parse >>> from psyclone.psyGen import PSyFactory >>> api = "gocean" >>> filename = "nemolite2d_alg.f90" >>> ast, invokeInfo = parse(filename, api=api, invoke_name="invoke") >>> psy = PSyFactory(api).create(invokeInfo) >>> >>> from psyclone.transformations import OMPParallelTrans, OMPSingleTrans >>> from psyclone.transformations import OMPTaskloopTrans >>> from psyclone.psyir.transformations import OMPTaskwaitTrans >>> singletrans = OMPSingleTrans() >>> paralleltrans = OMPParallelTrans() >>> tasklooptrans = OMPTaskloopTrans() >>> taskwaittrans = OMPTaskwaitTrans() >>> >>> schedule = psy.invokes.get('invoke_0').schedule >>> print(schedule.view()) >>> >>> # Apply the OpenMP Taskloop transformation to *every* loop >>> # in the schedule. >>> # This ignores loop dependencies. These are handled by the >>> # taskwait transformation. >>> for child in schedule.children: >>> tasklooptrans.apply(child, nogroup = true) >>> # 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) >>> taskwaittrans.apply(schedule.children) >>> print(schedule.view())
Inheritance

- apply(node, options=None)[source]#
Apply an OMPTaskwait Transformation to the supplied node (which must be an OMPParallelDirective). In the generated code this corresponds to adding zero or more OMPTaskwaitDirectives as appropriate:
!$OMP PARALLEL ... !$OMP TASKWAIT ... !$OMP TASKWAIT ... !$OMP END PARALLEL
- Parameters:
node (
psyclone.psyir.nodes.OMPParallelDirective) – the node to which to apply the transformation.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations and validation.
options["fail_on_no_taskloop"] (bool) – indicating whether this should throw an error if no OMPTaskloop nodes are found in this tree. This can be safely disabled as if there are no Taskloop nodes the result of this transformation is valid OpenMP code. Default is True
- static get_forward_dependence(taskloop, root)[source]#
Returns the next forward dependence for a taskloop using the dependence-analysis functionality provided by psyclone.psyir.tools.dependency_tools. Forward dependencies can be of the following types: Loop OMPDoDirective OMPTaskloopDirective OMPTaskwaitDirective (If in same OMPSingle and that single has nowait=False) OMPSingleDirective (If ancestor OMPSingle has nowait=False)
Loop, OMPDoDirective, OMPTaskloopDirective types are returned when a following directive is found which has a RaW, WaR or WaW dependency to taskloop.
An OMPTaskwaitDirective type is returned when a following directive is found inside the same parent OMPSingleDirective which has no nowait clause applied.
An OMPSingleDirective type is returned when the first dependency is within a different OMPSerialDirective, and the ancestor of taskloop is an OMPSingleDirective with no nowait clause.
The forward dependency is never a child of taskloop, and must have abs_position > taskloop.abs_position
- Parameters:
taskloop (
psyclone.psyir.nodes.OMPTaskloopDirective) – the taskloop node for which to find the forward_dependence.root (
psyclone.psyir.nodes.OMPParallelDirective) – the tree in which to search for the forward_dependence.
- Returns:
the forward_dependence of taskloop.
- Return type:
- validate(node, options=None)[source]#
Validity checks for input arguments.
- Parameters:
node (
psyclone.psyir.nodes.OMPParallelDirective) – the OMPParallelDirective node to validate.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
options["fail_on_no_taskloop"] (bool) – indicating whether this should throw an error if no OMPTaskloop nodes are found in this tree. This can be safely disabled as if there are no Taskloop nodes the result of this transformation is valid OpenMP code. Default is True.
- Raises:
TransformationError – If the supplied node is not an OMPParallelDirective
TransformationError – If there are no OMPTaskloopDirective nodes in this tree.
TransformationError – If taskloop dependencies can’t be satisfied due to dependencies across barrierless OpenMP Serial Regions.
- class psyclone.psyir.transformations.ParallelLoopTrans[source]#
Adds an abstract directive (it needs to be specified by sub-classing this transformation) to a loop indicating that it should be parallelised. It performs some data dependency checks to guarantee that the loop can be parallelised without changing the semantics of it.
If the nowait option is supplied to the apply function, then PSyclone will attempt to add the relevant asynchronous option to the directive, if supported. If its not supported, or PSyclone’s analysis suggests that it cannot be launched asynchronously, the transformation will occur as though the nowait option was not supplied. If the asynchronous option is added, PSyclone will also generate a corresponding barrier before the next dependent statement.
Inheritance

- apply(node, options=None, verbose=False, collapse=False, force=False, ignore_dependencies_for=None, privatise_arrays=False, sequential=False, nowait=False, reduction_ops=None, **kwargs)[source]#
Apply the Loop transformation to the specified node in a Schedule. This node must be a Loop since this transformation corresponds to wrapping the generated code with directives, e.g. for OpenMP:
!$OMP DO do ... ... end do !$OMP END DO
At code-generation time (when lowering is called), this node must be within (i.e. a child of) a PARALLEL region.
- Parameters:
node (
psyclone.psyir.nodes.Node) – the supplied node to which we will apply the loop parallelisation transformation.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
collapse (
Union[int,bool]) – if it’s a bool and is False (default), it won’t collapse. If it’s a bool and is True, it will collapse as much as possible. If it’s an integer, it will attempt to collapse until the specified number of loops (if they exist and are safe to collapse them). The options ‘ignore_dependencies_for’ and ‘force’ also affect the collapse applicability analysis.force (
bool) – whether to force parallelisation of the target loop (i.e. ignore any dependence analysis).ignore_dependencies_for (
Optional[List[str]]) – whether to ignore some symbol names from the dependence analysis checks.sequential (
bool) – whether this is a sequential loop.verbose (
bool) – whether to report the reasons the validate and collapse steps have failed.nowait (
bool) – whether to add a nowait clause and a corresponding barrier (or equivalent) to enable asynchronous execution.privatise_arrays (
bool) – whether to make the write after write dependency symbols declared as private.reduction_ops (
List[Union[Operator,Intrinsic]]) – if non-empty, attempt parallelisation of loops by inferring reduction clauses involving any of the reduction operators in the list.node_type_check (bool) – If the type of nodes enclosed in the loop should be tested to avoid including unsupported nodes in the transformation.
- validate(node, options=None, **kwargs)[source]#
Perform validation checks before applying the transformation
- Parameters:
node (
psyclone.psyir.nodes.Node) – the supplied node to which we will apply the loop parallelisation transformation.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
collapse (int | bool) – if it’s a bool and is False (default), it won’t collapse. If it’s a bool and is True, it will collapse as much as possible. If it’s an integer, it will attempt to collapse until the specified number of loops (if they exist and are safe to collapse them). The options ‘ignore_dependencies_for’ and ‘force’ also affect the collapse applicability analysis.
force (bool) – whether to force parallelisation of the target loop (i.e. ignore any dependence analysis).
ignore_dependencies_for (None | List[str]) – whether to ignore some symbol names from the dependence analysis checks.
sequential (bool) – whether this is a sequential loop.
verbose (bool) – whether to report the reasons the validate and collapse steps have failed.
nowait (bool) – whether to add a nowait clause and a corresponding barrier (or equivalent) to enable asynchronous execution.
privatise_arrays (bool) – whether to make the write after write dependency symbols declared as private.
reduction_ops (List[psyclone.psyir.nodes.operation.Operator | psyclone.psyir.nodes.intrinsic_call.IntrinsicCall.Intrinsic]) – if non-empty, attempt parallelisation of loops by inferring reduction clauses involving any of the reduction operators in the list.
node_type_check (bool) – If the type of nodes enclosed in the loop should be tested to avoid including unsupported nodes in the transformation.
reprod (bool) – indicating whether reproducible reductions should be used. By default the value from the config file will be used.
enable_reductions (bool) – whether to attempt to infer reduction clauses or not.
- Raises:
TypeError – if ‘collapse’ is not an int or a bool.
TypeError – if ‘ignore_dependencies_for’ is not a list of str.
TransformationError – if the given loop iterates over a colours (LFRic domain) iteration space.
TransformationError – if the given loops calls a procedure that is not guaranteed to be pure (and therefore could have dependencies beyond the specified by the arguments intent)
TransformationError – if the given loop is inside a pure routine as these do not allow parallel constructs.
TransformationError – if there is a data dependency that prevents the parallelisation of the loop and the provided options don’t disregard them.
- class psyclone.psyir.transformations.Product2LoopTrans[source]#
Provides a transformation from a PSyIR PRODUCT IntrinsicCall node to an equivalent PSyIR loop structure that is suitable for running in parallel on CPUs and GPUs. Validity checks are also performed.
If PRODUCT contains a single positional argument which is an array, the maximum value of all of the elements in the array is returned in the the scalar R.
R = PRODUCT(ARRAY)
For example, if the array is two dimensional, the equivalent code for real data is:
R = 1.0 DO J=LBOUND(ARRAY,2),UBOUND(ARRAY,2) DO I=LBOUND(ARRAY,1),UBOUND(ARRAY,1) R = R * ARRAY(I,J)
If the mask argument is provided then the mask is used to determine whether the product is applied:
R = PRODUCT(ARRAY, mask=MOD(ARRAY, 2.0)==1)
If the array is two dimensional, the equivalent code for real data is:
R = 1.0 DO J=LBOUND(ARRAY,2),UBOUND(ARRAY,2) DO I=LBOUND(ARRAY,1),UBOUND(ARRAY,1) IF (MOD(ARRAY(I,J), 2.0)==1) THEN R = R * ARRAY(I,J)
The dimension argument is currently not supported and will result in a TransformationError exception being raised.
R = PRODUCT(ARRAY, dimension=2)
The array passed to PRODUCT may use any combination of array syntax, array notation, array sections and scalar bounds:
R = PRODUCT(ARRAY) ! array syntax R = PRODUCT(ARRAY(:,:)) ! array notation R = PRODUCT(ARRAY(1:10,lo:hi)) ! array sections R = PRODUCT(ARRAY(1:10,:)) ! mix of array section and array notation R = PRODUCT(ARRAY(1:10,2)) ! mix of array section and scalar bound
An example use of this transformation is given below:
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.transformations import Product2LoopTrans >>> code = ("subroutine product_test(array)\n" ... " real :: array(10,10)\n" ... " real :: result\n" ... " result = product(array)\n" ... "end subroutine\n") >>> psyir = FortranReader().psyir_from_source(code) >>> product_node = psyir.children[0].children[0].children[1] >>> Product2LoopTrans().apply(product_node) >>> print(FortranWriter()(psyir)) subroutine product_test(array) real, dimension(10,10) :: array real :: result integer :: idx integer :: idx_1 result = 1.0 do idx = 1, 10, 1 do idx_1 = 1, 10, 1 result = result * array(idx_1,idx) enddo enddo end subroutine product_test
Inheritance

- apply(node, options=None, **kwargs)[source]#
Apply the Product2LoopTrans to the input node.
- Parameters:
node (
psyclone.psyir.nodes.IntrinsicCall) – the PRODUCT intrinsic to transform.options (Optional[Dict[str, Any]]) – options for the transformation.
- class psyclone.psyir.transformations.ProfileTrans[source]#
Create a profile region around a list of statements. For example:
>>> from psyclone.parse.algorithm import parse >>> from psyclone.parse.utils import ParseError >>> from psyclone.psyGen import PSyFactory, GenerationError >>> from psyclone.psyir.transformations import ProfileTrans >>> api = "gocean" >>> filename = "nemolite2d_alg.f90" >>> ast, invokeInfo = parse(filename, api=api, invoke_name="invoke") >>> psy = PSyFactory(api).create(invokeInfo) >>> >>> p_trans = ProfileTrans() >>> >>> schedule = psy.invokes.get('invoke_0').schedule >>> print(schedule.view()) >>> >>> # Enclose all children within a single profile region >>> p_trans.apply(schedule.children) >>> print(schedule.view())
This implementation relies completely on the base class PSyDataTrans for the actual work, it only adjusts the name etc, and the list of valid nodes.
Inheritance

- validate(nodes, options=None)[source]#
Checks that the supplied list of nodes is valid for profiling callipers.
- Parameters:
nodes (
psyclone.psyir.nodes.Nodeor list[psyclone.psyir.nodes.Node]) – a node or list of nodes to be instrumented with profiling.options["force"] (bool) – whether to ignore potential control flow jumps when applying this transformation. Default is False.
- Raises:
TransformationError – if the supplied region contains a potential control flow jump that could result in skipping the end of profiling caliper, e.g. EXIT or GOTO.
- class psyclone.psyir.transformations.PSyDataTrans(node_class=<class 'psyclone.psyir.nodes.psy_data_node.PSyDataNode'>)[source]#
Create a PSyData region around a list of statements. For example:
>>> from psyclone.parse.algorithm import parse >>> from psyclone.parse.utils import ParseError >>> from psyclone.psyGen import PSyFactory >>> api = "gocean" >>> ast, invoke_info = parse(SOURCE_FILE, api=api) >>> psy = PSyFactory(api).create(invoke_info) >>> >>> from psyclone.psyir.transformations import PSyDataTrans >>> data_trans = PSyDataTrans() >>> >>> schedule = psy.invokes.get('invoke_0').schedule >>> # Uncomment the following line to see a text view of the schedule >>> # print(schedule.view()) >>> >>> # Enclose all children within a single PSyData region >>> data_trans.apply(schedule.children) >>> # Uncomment the following line to see a text view of the schedule >>> # print(schedule.view()) >>> # Or to use custom region name: >>> data_trans.apply(schedule.children, ... {"region_name": ("module","region")})
- Parameters:
node_class (
psyclone.psyir.nodes.ExtractNode) – The Node class of which an instance will be inserted into the tree (defaults to PSyDataNode).
Inheritance

- apply(nodes, options=None)[source]#
Apply this transformation to a subset of the nodes within a schedule - i.e. enclose the specified Nodes in the schedule within a single PSyData region.
Note that if the nodes are within a routine that previously had the pure attribute, this attribute is removed.
- Parameters:
nodes (
psyclone.psyir.nodes.Nodeor list ofpsyclone.psyir.nodes.Node) – can be a single node or a list of nodes.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
options["prefix"] (str) – a prefix to use for the PSyData module name (
PREFIX_psy_data_mod) and the PSyDataType (PREFIX_PSYDATATYPE) - a “_” will be added automatically. It defaults to “”.options["region_name"] ((str,str)) – an optional name to use for this PSyData area, provided as a 2-tuple containing a location name followed by a local name. The pair of strings should uniquely identify a region unless aggregate information is required (and is supported by the runtime library).
- get_unique_region_name(nodes, options)[source]#
This function returns the region and module name. If they are specified in the user options, these names will just be returned (it is then up to the user to guarantee uniqueness). Otherwise a name based on the module and invoke will be created using indices to make sure the name is unique.
- Parameters:
nodes (list of
psyclone.psyir.nodes.Node) – a list of nodes.options (Dict[str, Any]) – a dictionary with options for transformations.
options["region_name"] ((str,str)) – an optional name to use for this PSyData area, provided as a 2-tuple containing a location name followed by a local name. The pair of strings should uniquely identify a region unless aggregate information is required (and is supported by the runtime library).
- property name#
This function returns the name of the transformation. It uses the Python 2/3 compatible way of returning the class name as a string, which means that the same function can be used for all derived classes.
- Returns:
the name of this transformation as a string.
- Return type:
str
- validate(nodes, options=None)[source]#
Checks that the supplied list of nodes is valid, that the location for this node is valid (not between a loop-directive and its loop), that there aren’t any name clashes with symbols that must be imported from the appropriate PSyData library and finally, calls the validate method of the base class.
- Parameters:
nodes ((list of)
psyclone.psyir.nodes.Loop) – a node or list of nodes to be instrumented with PSyData API calls.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
options["prefix"] (str) – a prefix to use for the PSyData module name (
PREFIX_psy_data_mod) and the PSyDataType (PREFIX_PSYDATATYPE) - a “_” will be added automatically. It defaults to “”.options["region_name"] ((str,str)) – an optional name to use for this PSyData area, provided as a 2-tuple containing a location name followed by a local name. The pair of strings should uniquely identify a region unless aggregate information is required (and is supported by the runtime library).
- Raises:
TransformationError – if the supplied list of nodes is empty.
TransformationError – if the PSyData node is inserted between an OpenMP/ACC directive and the loop(s) to which it applies.
TransformationError – if the ‘prefix’ or ‘region_name’ options are not valid.
TransformationError – if there will be a name clash between any existing symbols and those that must be imported from the appropriate PSyData library.
TransformationError – if the target nodes are within an ELEMENTAL routine.
- class psyclone.psyir.transformations.ReadOnlyVerifyTrans(node_class=<class 'psyclone.psyir.nodes.read_only_verify_node.ReadOnlyVerifyNode'>)[source]#
This transformation inserts a ReadOnlyVerifyNode or a node derived from ReadOnlyVerifyNode into the PSyIR of a schedule. At code creation time this node will use the PSyData API to create code that will verify that read-only quantities are not modified.
After applying the transformation the Nodes marked for verification are children of the ReadOnlyVerifyNode. Nodes to verify can be individual constructs within an Invoke (e.g. Loops containing a Kernel or BuiltIn call) or entire Invokes.
- Parameters:
node_class (
psyclone.psyir.nodes.ReadOnlyVerifyNodeor derived class) – The class of Node which will be inserted into the tree (defaults to ReadOnlyVerifyNode), but can be any derived class.
Inheritance

- validate(node_list, options=None)[source]#
Performs validation checks specific to read-only-based transformations.
- Parameters:
node_list (list of
psyclone.psyir.nodes.Node) – the list of Node(s) we are checking.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- Raises:
TransformationError – if transformation is applied to a Kernel or a BuiltIn call without its parent Loop.
TransformationError – if transformation is applied to a Loop without its parent Directive when optimisations are applied.
TransformationError – if transformation is applied to an orphaned Directive without its parent Directive.
- class psyclone.psyir.transformations.Reference2ArrayRangeTrans[source]#
Provides a transformation from PSyIR Array Notation (a reference to an Array) to a PSyIR Range. For example:
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Reference >>> from psyclone.psyir.transformations import TransformationError >>> CODE = ("program example\n" ... "real :: a(:)\n" ... "a = 0.0\n" ... "end program\n") >>> trans = Reference2ArrayRangeTrans() >>> psyir = FortranReader().psyir_from_source(CODE) >>> for reference in psyir.walk(Reference): ... try: ... trans.apply(reference) ... except TransformationError: ... pass >>> print(FortranWriter()(psyir)) program example real, dimension(:) :: a a(:) = 0.0 end program example
This transformation does not currently support arrays within structures, see issue #1858.
Inheritance

- apply(node, allow_call_arguments=False, **kwargs)[source]#
Apply the Reference2ArrayRangeTrans transformation to the specified node. The node must be a Reference to an array. The Reference is replaced by an ArrayReference with appropriate explicit range nodes (termed colon notation in Fortran).
- Parameters:
node (
psyclone.psyir.nodes.Reference) – a Reference node.allow_call_arguments (
bool) – by default, any references that may be arguments to non-elemental routines are not transformed. However, this transformation is sometimes used in other transformations where this restriction does not apply.
- validate(node, **kwargs)[source]#
Check that the node is a Reference node and that the symbol it references is an array.
- Parameters:
node (
psyclone.psyir.nodes.Reference) – a Reference node.allow_call_arguments (None) – by default, any references that may be arguments to non-elemental routines are not transformed. However, this transformation is sometimes used in other transformations where this restriction does not apply.
- Raises:
TransformationError – if the node is not a Reference node or the Reference node not does not reference an array symbol.
TransformationError – if the Reference node is (or may be) passed as an argument to a call that is not elemental and allow_call_arguments is False.
- class psyclone.psyir.transformations.RegionTrans[source]#
This abstract class is a base class for all transformations that act on a list of nodes. It gives access to a validate function that makes sure that the nodes in the list are in the same order as in the original AST, no node is duplicated, and that all nodes have the same parent. We also check that all nodes to be enclosed are valid for this transformation - this requires that the sub-class populate the excluded_node_types tuple.
Inheritance

- apply(nodes, node_type_check=True, options=None, **kwargs)[source]#
Apply the RegionTrans transformation.
- Parameters:
nodes (
Union[Node,Schedule,List[Node]]) – can be a single node, a schedule or a list of nodes.node_type_check (
bool) – whether or not the type of the nodes enclosed in the region should be tested to avoid using unsupported nodes inside a region.node (depends on actual transformation) – The node (or list of nodes) for the transformation - specific to the actual transform used.
- get_node_list(nodes)[source]#
This is a helper function for region based transformations. The parameter for any of those transformations is either a single Node, a Schedule, or a list of nodes. This function converts this into a list of nodes according to the parameter type. This function will always return a copy, to avoid issues e.g. if a child list of a node should be provided, and a transformation changes the order in this list (which would then also change the order of the nodes in the tree).
If nodes happens to be a list containing a single Schedule node then the behaviour is the same as if it had been a single Schedule node, i.e. we return a list of that Schedule’s children.
- Parameters:
nodes (
Union[Node,Schedule,List[Node]]) – can be a single node, a schedule or a list of nodes.- Return type:
List[Node]- Returns:
a list of nodes.
- Raises:
TransformationError – if the supplied parameter is neither a single Node, nor a Schedule, nor a list of Nodes.
- validate(nodes, options=None, **kwargs)[source]#
Checks that the nodes in nodes are valid for a region transformation.
- Parameters:
nodes (
Union[Node,Schedule,List[Node]]) – can be a single node, a schedule or a list of nodes.node_type_check (bool) – whether or not the type of the nodes enclosed in the region should be tested to avoid using unsupported nodes inside a region.
node (depends on actual transformation) – The node (or list of nodes) for the transformation - specific to the actual transform used.
- Raises:
TransformationError – if the nodes in the list are not in the original order in which they are in the AST, a node is duplicated or the nodes have different parents.
TransformationError – if any of the nodes to be enclosed in the region are of an unsupported type.
TransformationError – if the parent of the supplied Nodes is not a Schedule or a Directive.
TransformationError – if the supplied options are not a dictionary.
- class psyclone.psyir.transformations.ReplaceInductionVariablesTrans[source]#
Move all supported induction variables out of the loop, and replace their usage inside the loop. For example:
>>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Loop >>> from psyclone.psyir.transformations import ReplaceInductionVariablesTrans >>> from psyclone.psyir.backend.fortran import FortranWriter >>> psyir = FortranReader().psyir_from_source(""" ... subroutine sub() ... integer :: i, im, ic, tmp(100) ... do i=1, 100 ... im = i - 1 ... ic = 2 ... tmp(i) = ic * im ... enddo ... end subroutine sub""") >>> loop = psyir.walk(Loop)[0] >>> ReplaceInductionVariablesTrans().apply(loop) >>> print(FortranWriter()(psyir)) subroutine sub() integer :: i integer :: im integer :: ic integer, dimension(100) :: tmp do i = 1, 100, 1 tmp(i) = 2 * (i - 1) enddo ic = 2 im = i - 1 - 1 end subroutine sub
The replaced induction variables assignments are added after the loop, so these variables will have the correct value if they are used elsewhere.
The following restrictions apply for the assignment to an induction variable:
the variable must be a scalar (i.e. no array access at all, not even a constant like a(3) or a%b(3)%c)
none of variables on the right-hand side can be written in the loop body (the loop variable is written in the Loop statement, not in the body, so it can be used).
Only intrinsic function calls are allowed on the RHS (since they are known to be elemental)
the assigned variable must not be read before the assignment.
the assigned variable cannot occur on the right-hand side (e.g. k = k + 3).
there must be only one assignment to this induction variable.
Inheritance

- apply(node, options=None)[source]#
Apply the ReplaceInductionVariablesTrans transformation to the specified node. The node must be a loop. In case of nested loops, the transformation might need to be applied several times, from the inner-most loop outwards.
- Parameters:
node (
psyclone.psyir.nodes.Loop) – a Loop node.
- validate(node, options=None)[source]#
Perform various checks to ensure that it is valid to apply the ReplaceInductionVariablesTrans transformation to the supplied PSyIR Node.
- Parameters:
node (
psyclone.psyir.nodes.Assignment) – the node that is being checked.- Raises:
TransformationError – if the node argument is not a Loop.
- class psyclone.psyir.transformations.ReplaceReferenceByLiteralTrans[source]#
This transformation takes a psyir Routine and replace all Reference psyir Nodes by Literal if the corresponding symbol from the symbol table is constant. That is to say the symbol is a Fortran parameter. NOTE: this transformation does not handle parent scopes recursively yet. For example:
>>> from psyclone.psyir.backend.fortran import FortranWriter >>> from psyclone.psyir.symbols import INTEGER_TYPE >>> from psyclone.psyir.transformations import ( ReplaceReferenceByLiteralTrans) >>> source = """program test ... use mymod ... type(my_type):: t1, t2, t3, t4 ... integer, parameter :: x=3, y=12, z=13 ... integer, parameter :: u1=1, u2=2, u3=3, u4=4 ... integer i, invariant, ic1, ic2, ic3 ... real, dimension(10) :: a ... invariant = 1 ... do i = 1, 10 ... t1%a = z ... a(ic1) = u1+(ic1+x)*ic1 ... a(ic2) = u2+(ic2+y)*ic2 ... a(ic3) = u3+(ic3+z)*ic3 ... a(t1%a) = u4+(t1%a+u4*z)*t1%a ... end do ... end program test""" >>> fortran_writer = FortranWriter() >>> fortran_reader = FortranReader() >>> psyir = fortran_reader.psyir_from_source(source) >>> routine = psyir.walk(Routine)[0] >>> rrbl = ReplaceReferenceByLiteralTrans() >>> rrbl.apply(routine) >>> written_code = fortran_writer(routine) >>> print(written_code) program test use mymod integer, parameter :: x = 3 integer, parameter :: y = 12 integer, parameter :: z = 13 integer, parameter :: u1 = 1 integer, parameter :: u2 = 2 integer, parameter :: u3 = 3 integer, parameter :: u4 = 4 type(my_type) :: t1 type(my_type) :: t2 type(my_type) :: t3 type(my_type) :: t4 integer :: i integer :: invariant integer :: ic1 integer :: ic2 integer :: ic3 real, dimension(10) :: a invariant = 1 do i = 1, 10, 1 t1%a = 13 a(ic1) = 1 + (ic1 + 3) * ic1 a(ic2) = 2 + (ic2 + 12) * ic2 a(ic3) = 3 + (ic3 + 13) * ic3 a(t1%a) = 4 + (t1%a + 4 * 13) * t1%a enddo end program test
Inheritance

- apply(node, options=None)[source]#
Applies the transformation to a Routine node: * First update a dictionary (param_table) with the Literal of constant (parameter) symbol from node.parent symbol_table, and from node.symbol_table. * Second, use this updated param_table to replace reference in node psyir_tree with the corresponsing Literal. * Third, use this updated param_table to replace reference in node symbol_table DataSymbol array’s dimensions with the corresponding Literal.
Note: If an update cannot be performed for any reason then the substitution is skipped and a comment starting with ‘Psyclone(ReplaceReferenceByLiteral):’ is added to the generated code.
- Parameters:
node (
Routine) – node on which the transformation is appliedoptions (
Optional[Dict[str,Any]]) – a dictionary with options for transformations.
- validate(node, options=None)[source]#
Perform various checks to ensure that it is valid to apply the ReplaceReferenceByLiteralTrans transformation to the supplied PSyIR Node.
- Parameters:
node – the node that is being checked.
options (_type_, optional) – not used, defaults to None
- Raises:
TransformationError – if the node argument is not a Routine.
- class psyclone.psyir.transformations.ValueRangeCheckTrans(node_class=<class 'psyclone.psyir.nodes.value_range_check_node.ValueRangeCheckNode'>)[source]#
This transformation inserts a ValueRangeCheckNode into the PSyIR of a schedule. At code creation time this node will use the PSyData API to create code that will verify all input parameters are not NANs and not infinite, and the same for all output parameters.
After applying the transformation the Nodes marked for verification are children of the ValueRangeCheckNode. Nodes to verify can be individual constructs within an Invoke (e.g. Loops containing a Kernel or BuiltIn call) or entire Invokes.
- Parameters:
node_class (
psyclone.psyir.nodes.ValueRangeCheckNodeor derived class) – The class of Node which will be inserted into the tree (defaults to ValueRangeCheckNode), but can be any derived class.
Inheritance

- validate(node_list, options=None)[source]#
Performs validation checks specific to nan-test transformations. This function is only here so that it is documented.
- Parameters:
node_list (list of
psyclone.psyir.nodes.Node) – the list of Node(s) we are checking.options (Optional[Dict[str, Any]]) – a dictionary with options for transformations.
- Raises:
TransformationError – if transformation is applied to a Kernel or a BuiltIn call without its parent Loop.
TransformationError – if transformation is applied to a Loop without its parent Directive when optimisations are applied.
TransformationError – if transformation is applied to an orphaned Directive without its parent Directive.
- class psyclone.psyir.transformations.ParallelRegionTrans[source]#
Base class for transformations that create a parallel region.
Inheritance

- 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 nodes are not all children of the same parent (siblings).
- class psyclone.psyir.transformations.OMPTaskloopTrans(grainsize=None, num_tasks=None, nogroup=False)[source]#
Adds an OpenMP taskloop directive to a loop. Only one of grainsize or num_tasks must be specified.
TODO: #1364 Taskloops do not yet support reduction clauses.
- Parameters:
grainsize (int or None) – the grainsize to use in for this transformation.
num_tasks (int or None) – the num_tasks to use for this transformation.
nogroup (bool) – whether or not to use a nogroup clause for this transformation. Default is False.
For example:
>>> from pysclone.parse.algorithm import parse >>> from psyclone.psyGen import PSyFactory >>> api = "gocean" >>> ast, invokeInfo = parse(GOCEAN_SOURCE_FILE, api=api) >>> psy = PSyFactory(api).create(invokeInfo) >>> >>> from psyclone.transformations import OMPParallelTrans, OMPSingleTrans >>> from psyclone.transformations import OMPTaskloopTrans >>> from psyclone.psyir.transformations import OMPTaskwaitTrans >>> singletrans = OMPSingleTrans() >>> paralleltrans = OMPParallelTrans() >>> tasklooptrans = OMPTaskloopTrans() >>> taskwaittrans = OMPTaskwaitTrans() >>> >>> 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 Taskloop transformation to *every* loop >>> # in the schedule. >>> # This ignores loop dependencies. These can be handled >>> # by the OMPTaskwaitTrans >>> for child in schedule.children: >>> tasklooptrans.apply(child) >>> # 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) >>> # Ensure loop dependencies are satisfied >>> taskwaittrans.apply(schedule.children) >>> # Uncomment the following line to see a text view of the schedule >>> # print(schedule.view())
Inheritance

- apply(node, options=None, **kwargs)[source]#
Apply the OMPTaskloopTrans transformation to the specified node in a Schedule. This node must be a Loop since this transformation corresponds to wrapping the generated code with directives like so:
!$OMP TASKLOOP do ... ... end do !$OMP END TASKLOOP
At code-generation time (when lowering is called), this node must be within (i.e. a child of) an OpenMP SERIAL region.
If the keyword “nogroup” is specified in the options, it will cause a nogroup clause be generated if it is set to True. This will override the value supplied to the constructor, but will only apply to the apply call to which the value is supplied.
- Parameters:
node (
psyclone.psyir.nodes.Node) – the supplied node to which we will apply the OMPTaskloopTrans transformationoptions (Optional[Dict[str, Any]]) – a dictionary with options for transformations and validation.
options["nogroup"] (bool) – indicating whether a nogroup clause should be applied to this taskloop.
- property omp_grainsize#
Returns the grainsize that will be specified by this transformation. By default the grainsize clause is not applied, so grainsize is None.
- Returns:
The grainsize specified by this transformation.
- Return type:
int or None
- property omp_nogroup#
Returns whether the nogroup clause should be specified for this transformation. By default the nogroup clause is applied.
- Returns:
whether the nogroup clause should be specified by this transformation.
- Return type:
bool
- property omp_num_tasks#
Returns the num_tasks that will be specified by this transformation. By default the num_tasks clause is not applied so num_tasks is None.
- Returns:
The grainsize specified by this transformation.
- Return type:
int or None
- class psyclone.psyir.transformations.OMPDeclareTargetTrans[source]#
Adds an OpenMP declare target directive to the specified routine.
For example:
>>> from psyclone.psyir.frontend.fortran import FortranReader >>> from psyclone.psyir.nodes import Loop >>> from psyclone.transformations import OMPDeclareTargetTrans >>> >>> tree = FortranReader().psyir_from_source(""" ... subroutine my_subroutine(A) ... integer, dimension(10, 10), intent(inout) :: A ... integer :: i ... integer :: j ... do i = 1, 10 ... do j = 1, 10 ... A(i, j) = 0 ... end do ... end do ... end subroutine ... """ >>> omptargettrans = OMPDeclareTargetTrans() >>> omptargettrans.apply(tree.walk(Routine)[0])
will generate:
subroutine my_subroutine(A) integer, dimension(10, 10), intent(inout) :: A integer :: i integer :: j !$omp declare target do i = 1, 10 do j = 1, 10 A(i, j) = 0 end do end do end subroutine
Inheritance

- apply(node, options=None)[source]#
Insert an OMPDeclareTargetDirective inside the provided routine or associated PSyKAl kernel.
- Parameters:
node (
psyclone.psyir.nodes.Routine|psyclone.psyGen.Kern) – 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.
options["device_string"] (str) – provide a compiler-platform identifier.
- validate(node, options=None)[source]#
Check that an OMPDeclareTargetDirective can be inserted.
- 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.
options["device_string"] (str) – provide a compiler-platform identifier.
- 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.
Exceptions#
TransformationError: Provides a PSyclone-specific error class for errors found during
