loki.transformations.block_index_transformations

Classes

BlockViewToFieldViewTransformation(horizontal)

A very IFS-specific transformation to replace per-block, i.e. per OpenMP-thread, view pointers with per-field view pointers.

InjectBlockIndexTransformation(block_dim)

A transformation pass to inject the block-index in arrays promoted by a previous transformation pass.

LowerBlockIndexTransformation(block_dim[, ...])

Transformation to lower the block index via appending the block index to variable dimensions/shape.

LowerBlockLoopTransformation(block_dim)

Lower the block loop to calls within this loop.

class BlockViewToFieldViewTransformation(horizontal, global_gfl_ptr=False)

Bases: Transformation

A very IFS-specific transformation to replace per-block, i.e. per OpenMP-thread, view pointers with per-field view pointers. It should be noted that this transformation only replaces the view pointers but does not actually insert the block index into the promoted view pointers. Therefore this transformation must always be followed by the InjectBlockIndexTransformation.

For example, the following code:

do jlon=1,nproma
  mystruct%p(jlon,:) = 0.
enddo

is transformed to:

do jlon=1,nproma
  mystruct%p_field(jlon,:) = 0.
enddo

As the rank of my_struct%p_field is one greater than that of my_struct%p, we would need to also apply the InjectBlockIndexTransformation to obtain semantically correct code:

do jlon=1,nproma
  mystruct%p_field(jlon,:,ibl) = 0.
enddo

Specific arrays in individual routines can also be marked for exclusion from this transformation by assigning them to the exclude_arrays list in the SchedulerConfig.

This transformation also creates minimal definitions of FIELD API wrappers (i.e. FIELD_RANKSUFF_ARRAY) and uses them to enrich the DataType of relevant variable declarations and expression nodes. This is required because FIELD API can be built independently of library targets Loki would typically operate on.

Parameters:
  • horizontal (Dimension) – Dimension object describing the variable conventions used in code to define the horizontal data dimension and iteration space.

  • global_gfl_ptr (bool) – Toggle whether thread-local gfl_ptr should be replaced with global.

  • key (str, optional) – Specify a different identifier under which trafo_data is stored

item_filter = (<class 'loki.batch.item.ProcedureItem'>,)
transform_subroutine(routine, **kwargs)

Defines the transformation to apply to Subroutine items.

For transformations that modify Subroutine objects, this method should be implemented. It gets called via the dispatch method apply().

Parameters:
  • routine (Subroutine) – The subroutine to be transformed.

  • **kwargs (optional) – Keyword arguments for the transformation.

static propagate_defs_to_children(key, definitions, successors)

Enrich all successors with the dummy FIELD_API definitions.

process_driver(routine, successors)
build_ydvars_global_gfl_ptr(var)

Replace accesses to thread-local YDVARS%GFL_PTR with global YDVARS%GFL_PTR_G.

process_body(body, definitions, successors, targets, exclude_arrays)
process_kernel(routine, item, successors, targets, exclude_arrays)
class InjectBlockIndexTransformation(block_dim)

Bases: Transformation

A transformation pass to inject the block-index in arrays promoted by a previous transformation pass. As such, this transformation also relies on the block-index, or a known alias, being already present in routines that are to be transformed.

For array access in a Subroutine body, it operates by comparing the local shape of an array with its declared shape. If the local shape is of rank one less than the declared shape, then the block-index is appended to the array’s dimensions.

For CallStatement arguments, if the rank of the argument is one less than that of the corresponding dummy-argument, the block-index is appended to the argument’s dimensions. It should be noted that this logic relies on the CallStatement being free of any sequence-association.

For example, the following code:

subroutine kernel1(nblks, ...)
   ...
   integer, intent(in) :: nblks
   integer :: ibl
   real :: var(jlon,nlev,nblks)

   do ibl=1,nblks
     do jlon=1,nproma
       var(jlon,:) = 0.
     enddo

     call kernel2(var,...)
   enddo
   ...
end subroutine kernel1

subroutine kernel2(var, ...)
   ...
   real :: var(jlon,nlev)
end subroutine kernel2

is transformed to:

subroutine kernel1(nblks, ...)
   ...
   integer, intent(in) :: nblks
   integer :: ibl
   real :: var(jlon,nlev,nblks)

   do ibl=1,nblks
     do jlon=1,nproma
       var(jlon,:,ibl) = 0.
     enddo

     call kernel2(var(:,:,ibl),...)
   enddo
   ...
end subroutine kernel1

subroutine kernel2(var, ...)
   ...
   real :: var(jlon,nlev)
end subroutine kernel2

Specific arrays in individual routines can also be marked for exclusion from this transformation by assigning them to the exclude_arrays list in the SchedulerConfig.

Parameters:
  • block_dim (Dimension) – Dimension object describing the variable conventions used in code to define the blocking data dimension and iteration space.

  • key (str, optional) – Specify a different identifier under which trafo_data is stored

item_filter = (<class 'loki.batch.item.ProcedureItem'>,)
transform_subroutine(routine, **kwargs)

Defines the transformation to apply to Subroutine items.

For transformations that modify Subroutine objects, this method should be implemented. It gets called via the dispatch method apply().

Parameters:
  • routine (Subroutine) – The subroutine to be transformed.

  • **kwargs (optional) – Keyword arguments for the transformation.

static get_call_arg_rank(arg)

Utility to retrieve the local rank of a CallStatement argument.

get_block_index(routine)

Utility to retrieve the block-index loop induction variable.

process_body(body, block_index, targets, exclude_arrays)
process_kernel(routine, targets, exclude_arrays)
class LowerBlockIndexTransformation(block_dim, recurse_to_kernels=True)

Bases: Transformation

Transformation to lower the block index via appending the block index to variable dimensions/shape. However, this only handles variable declarations/definitions. Therefore this transformation must always be followed by the InjectBlockIndexTransformation.

For example, the following code:

SUBROUTINE driver (nlon, nlev, nb, var)
    INTEGER, INTENT(IN) :: nlon, nlev, nb
    REAL, INTENT(INOUT) :: var(nlon, nlev, nb)
    DO ibl=1,10
        CALL kernel(nlon, nlev, var(:, :, ibl))
    END DO
END SUBROUTINE driver

SUBROUTINE kernel (nlon, nlev, var, another_var, icend, lstart, lend)
    INTEGER, INTENT(IN) :: nlon, nlev
    REAL, INTENT(INOUT) :: var(nlon, nlev)
    DO jk=1,nlev
        DO jl=1,nlon
            var(jl, jk) = 0.
        END DO
    END DO
END SUBROUTINE kernel

is transformed to:

SUBROUTINE driver (nlon, nlev, nb, var)
    INTEGER, INTENT(IN) :: nlon, nlev, nb
    REAL, INTENT(INOUT) :: var(nlon, nlev, nb)
    DO ibl=1,10
        CALL kernel(nlon, nlev, var(:, :, :), ibl=ibl, nb=nb)
    END DO
END SUBROUTINE driver

SUBROUTINE kernel (nlon, nlev, var, another_var, icend, lstart, lend, ibl, nb)
    INTEGER, INTENT(IN) :: nlon, nlev, ibl, nb
    REAL, INTENT(INOUT) :: var(nlon, nlev, nb)
    DO jk=1,nlev
        DO jl=1,nlon
            var(jl, jk) = 0.
        END DO
    END DO
END SUBROUTINE kernel

Warning

The block index is not injected by this transformation. To inject the block index to the promoted arrays call the InjectBlockIndexTransformation after this transformation!

Parameters:
  • block_dim (Dimension) – Dimension object describing the variable conventions used in code to define the blocking data dimension and iteration space.

  • recurse_to_kernels (bool, optional) – Recurse/continue with/to (nested) kernels and lower the block index for those as well (default: False).

item_filter = (<class 'loki.batch.item.ProcedureItem'>,)
transform_subroutine(routine, **kwargs)

Defines the transformation to apply to Subroutine items.

For transformations that modify Subroutine objects, this method should be implemented. It gets called via the dispatch method apply().

Parameters:
  • routine (Subroutine) – The subroutine to be transformed.

  • **kwargs (optional) – Keyword arguments for the transformation.

process(routine, targets, role)

This method adds the blocking index and size as arguments (if not already there yet) for calls and corresponding callees and updates the dimension and shape of arguments for calls and callees.

Note

The injection of the block index for promoted arrays is not part of this method (and also not part of this transformation!)

class LowerBlockLoopTransformation(block_dim)

Bases: Transformation

Lower the block loop to calls within this loop.

For example, the following code:

subroutine driver(nblks, ...)
    ...
    integer, intent(in) :: nblks
    integer :: ibl
    real :: var(jlon,nlev,nblks)

    do ibl=1,nblks
        call kernel2(var,...,nblks,ibl)
    enddo
    ...
end subroutine driver

subroutine kernel(var, ..., nblks, ibl)
    ...
    real :: var(jlon,nlev,nblks)

    do jl=1,...
        do jk=1,...
            var(jk,jl,ibl) = ...
        end do
    end do
end subroutine kernel

is transformed to:

subroutine driver(nblks, ...)
    ...
    integer, intent(in) :: nblks
    integer :: ibl
    real :: var(jlon,nlev,nblks)

    call kernel2(var,..., nblks)
    ...
end subroutine driver

subroutine kernel(var, ..., nblks)
    ...
    integer :: ibl
    real :: var(jlon,nlev,nblks)

    do ibl=1,nblks
        do jl=1,...
            do jk=1,...
                var(jk,jl,ibl) = ...
            end do
        end do
    end do
end subroutine kernel
Parameters:

block_dim (Dimension) – Dimension object describing the variable conventions used in code to define the blocking data dimension and iteration space.

item_filter = (<class 'loki.batch.item.ProcedureItem'>,)
transform_subroutine(routine, **kwargs)

Defines the transformation to apply to Subroutine items.

For transformations that modify Subroutine objects, this method should be implemented. It gets called via the dispatch method apply().

Parameters:
  • routine (Subroutine) – The subroutine to be transformed.

  • **kwargs (optional) – Keyword arguments for the transformation.

static arg_to_local_var(routine, var)
local_var(call, variables)
static generate_pragma(loop)
update_call_signature(call, loop, loop_defined_symbols, additional_kwargs)
process_driver(routine, targets)