Skip to content

Commit

Permalink
Extracts: Rename extract_marked_routines to outline_pragma_regions
Browse files Browse the repository at this point in the history
This also changes the semantics of the associated pragma annotations
to `!$loki outline`.
  • Loading branch information
mlange05 committed Oct 14, 2024
1 parent cc85694 commit 4663da6
Show file tree
Hide file tree
Showing 4 changed files with 91 additions and 92 deletions.
22 changes: 11 additions & 11 deletions loki/transformations/extract/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -17,7 +17,7 @@
"""

from loki.transformations.extract.internal import * # noqa
from loki.transformations.extract.marked import * # noqa
from loki.transformations.extract.outline import * # noqa

from loki.batch import Transformation

Expand All @@ -33,16 +33,16 @@ class ExtractTransformation(Transformation):
Parameters
----------
inline_internals : bool
extract_internals : bool
Extract internal procedure (see :any:`extract_internal_procedures`);
default: False.
inline_marked : bool
Extract :any:`Subroutine` objects marked by pragma annotations
(see :any:`extract_marked_subroutines`); default: True.
outline_regions : bool
Outline pragma-annotated code regions to :any:`Subroutine` objects.
(see :any:`outline_pragma_regions`); default: True.
"""
def __init__(self, extract_internals=False, extract_marked=True):
def __init__(self, extract_internals=False, outline_regions=True):
self.extract_internals = extract_internals
self.extract_marked = extract_marked
self.outline_regions = outline_regions

def transform_module(self, module, **kwargs):
"""
Expand All @@ -57,9 +57,9 @@ def transform_module(self, module, **kwargs):
module.contains.append(new_routines)

# Extract pragma-marked code regions into standalone subroutines
if self.extract_marked:
if self.outline_regions:
for routine in module.subroutines:
new_routines = extract_marked_subroutines(routine)
new_routines = outline_pragma_regions(routine)
module.contains.append(new_routines)

def transform_file(self, sourcefile, **kwargs):
Expand All @@ -75,7 +75,7 @@ def transform_file(self, sourcefile, **kwargs):
sourcefile.ir.append(new_routines)

# Extract pragma-marked code regions into standalone subroutines
if self.extract_marked:
if self.outline_regions:
for routine in sourcefile.subroutines:
new_routines = extract_marked_subroutines(routine)
new_routines = outline_pragma_regions(routine)
sourcefile.ir.append(new_routines)
Original file line number Diff line number Diff line change
Expand Up @@ -17,17 +17,16 @@
from loki.tools import as_tuple, CaseInsensitiveDict


__all__ = ['extract_marked_subroutines']
__all__ = ['outline_pragma_regions']


def extract_marked_subroutines(routine):
def outline_pragma_regions(routine):
"""
Convert regions annotated with ``!$loki extract`` pragmas to subroutine calls.
Convert regions annotated with ``!$loki outline`` pragmas to subroutine calls.
The pragma syntax for regions to convert to subroutines is
``!$loki extract [name(...)] [in(...)] [out(...)] [inout(...)]``
and ``!$loki end extract``.
``!$loki outline [name(...)] [in(...)] [out(...)] [inout(...)]``
and ``!$loki end outline``.
A new subroutine is created with the provided name (or an auto-generated default name
derived from the current subroutine name) and the content of the pragma region as body.
Expand All @@ -51,12 +50,12 @@ def extract_marked_subroutines(routine):
with pragma_regions_attached(routine):
with dataflow_analysis_attached(routine):
for region in FindNodes(PragmaRegion).visit(routine.body):
if not is_loki_pragma(region.pragma, starts_with='extract'):
if not is_loki_pragma(region.pragma, starts_with='outline'):
continue

# Name the external routine
parameters = get_pragma_parameters(region.pragma, starts_with='extract')
name = parameters.get('name', f'{routine.name}_extracted_{counter}')
parameters = get_pragma_parameters(region.pragma, starts_with='outline')
name = parameters.get('name', f'{routine.name}_outlined_{counter}')
counter += 1

# Create the external subroutine containing the routine's imports and the region's body
Expand Down
34 changes: 17 additions & 17 deletions loki/transformations/extract/tests/test_extract_transformation.py
Original file line number Diff line number Diff line change
Expand Up @@ -15,9 +15,9 @@


@pytest.mark.parametrize('frontend', available_frontends())
@pytest.mark.parametrize('extract_marked', [False, True])
@pytest.mark.parametrize('outline_regions', [False, True])
@pytest.mark.parametrize('extract_internals', [False, True])
def test_extract_transformation_module(extract_internals, extract_marked, frontend):
def test_extract_transformation_module(extract_internals, outline_regions, frontend):
"""
Test basic subroutine extraction from marker pragmas in modules.
"""
Expand All @@ -37,13 +37,13 @@ def test_extract_transformation_module(extract_internals, extract_marked, fronte
y(i,:) = b(i)
end do
!$loki extract name(test1)
!$loki outline name(test1)
do i=1, n
do j=1, n+1
x(i) = x(i) + y(i, j)
end do
end do
!$loki end extract
!$loki end outline
do i=1, n
call plus_one(x, i=i)
Expand All @@ -62,26 +62,26 @@ def test_extract_transformation_module(extract_internals, extract_marked, fronte
module = Module.from_source(fcode, frontend=frontend)

ExtractTransformation(
extract_internals=extract_internals, extract_marked=extract_marked
extract_internals=extract_internals, outline_regions=outline_regions
).apply(module)

routines = tuple(r for r in module.contains.body if isinstance(r, Subroutine))
assert len(routines) == 1 + (1 if extract_internals else 0) + (1 if extract_marked else 0)
assert len(routines) == 1 + (1 if extract_internals else 0) + (1 if outline_regions else 0)
assert ('plus_one' in module) == extract_internals
assert ('test1' in module) == extract_marked
assert ('test1' in module) == outline_regions

outer = module['outer']
assert len(FindNodes(ir.CallStatement).visit(outer.body)) == (2 if extract_marked else 1)
assert len(FindNodes(ir.CallStatement).visit(outer.body)) == (2 if outline_regions else 1)
outer_internals = tuple(r for r in outer.contains.body if isinstance(r, Subroutine))
assert len(outer_internals) == (0 if extract_internals else 1)


@pytest.mark.parametrize('frontend', available_frontends())
@pytest.mark.parametrize('extract_marked', [False, True])
@pytest.mark.parametrize('outline_regions', [False, True])
@pytest.mark.parametrize('extract_internals', [False, True])
def test_extract_transformation_sourcefile(extract_internals, extract_marked, frontend):
def test_extract_transformation_sourcefile(extract_internals, outline_regions, frontend):
"""
Test internal procedure extraction from subroutines.
Test internal procedure extraction and region outlining from subroutines.
"""
fcode = """
subroutine outer(n, a, b)
Expand All @@ -95,13 +95,13 @@ def test_extract_transformation_sourcefile(extract_internals, extract_marked, fr
y(i,:) = b(i)
end do
!$loki extract name(test1)
!$loki outline name(test1)
do i=1, n
do j=1, n+1
x(i) = x(i) + y(i, j)
end do
end do
!$loki end extract
!$loki end outline
do i=1, n
call plus_one(x, i=i)
Expand All @@ -119,15 +119,15 @@ def test_extract_transformation_sourcefile(extract_internals, extract_marked, fr
sourcefile = Sourcefile.from_source(fcode, frontend=frontend)

ExtractTransformation(
extract_internals=extract_internals, extract_marked=extract_marked
extract_internals=extract_internals, outline_regions=outline_regions
).apply(sourcefile)

routines = tuple(r for r in sourcefile.ir.body if isinstance(r, Subroutine))
assert len(routines) == 1 + (1 if extract_internals else 0) + (1 if extract_marked else 0)
assert len(routines) == 1 + (1 if extract_internals else 0) + (1 if outline_regions else 0)
assert ('plus_one' in sourcefile) == extract_internals
assert ('test1' in sourcefile) == extract_marked
assert ('test1' in sourcefile) == outline_regions

outer = sourcefile['outer']
assert len(FindNodes(ir.CallStatement).visit(outer.body)) == (2 if extract_marked else 1)
assert len(FindNodes(ir.CallStatement).visit(outer.body)) == (2 if outline_regions else 1)
outer_internals = tuple(r for r in outer.contains.body if isinstance(r, Subroutine))
assert len(outer_internals) == (0 if extract_internals else 1)
Loading

0 comments on commit 4663da6

Please sign in to comment.