Reference Guide  2.5.0
kernel_module_inline_trans.py
1 # -----------------------------------------------------------------------------
2 # BSD 3-Clause License
3 #
4 # Copyright (c) 2017-2024, Science and Technology Facilities Council.
5 # All rights reserved.
6 #
7 # Redistribution and use in source and binary forms, with or without
8 # modification, are permitted provided that the following conditions are met:
9 #
10 # * Redistributions of source code must retain the above copyright notice, this
11 # list of conditions and the following disclaimer.
12 #
13 # * Redistributions in binary form must reproduce the above copyright notice,
14 # this list of conditions and the following disclaimer in the documentation
15 # and/or other materials provided with the distribution.
16 #
17 # * Neither the name of the copyright holder nor the names of its
18 # contributors may be used to endorse or promote products derived from
19 # this software without specific prior written permission.
20 #
21 # THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
22 # "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
23 # LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
24 # FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
25 # COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
26 # INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
27 # BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
28 # LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
29 # CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
30 # LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
31 # ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
32 # POSSIBILITY OF SUCH DAMAGE.
33 # -----------------------------------------------------------------------------
34 # Authors R. W. Ford, A. R. Porter, S. Siso and N. Nobre, STFC Daresbury Lab
35 # A. B. G. Chalk STFC Daresbury Lab
36 # J. Henrichs, Bureau of Meteorology
37 # Modified I. Kavcic, Met Office
38 
39 ''' This module provides the KernelModuleInlineTrans transformation. '''
40 
41 
42 from psyclone.psyGen import Transformation, CodedKern
43 from psyclone.psyir.transformations import TransformationError
44 from psyclone.psyir.symbols import (
45  RoutineSymbol, DataSymbol, IntrinsicSymbol, DataTypeSymbol, Symbol,
46  ContainerSymbol, DefaultModuleInterface)
47 from psyclone.psyir.nodes import (
48  Container, ScopingNode, Reference, Routine, Literal, CodeBlock, Call)
49 
50 
52  ''' Module-inlines (bring the subroutine to the same compiler-unit) the
53  subroutine pointed by this Kernel. For example:
54 
55  .. code-block:: python
56 
57  from psyclone.domain.common.transformations import \\
58  KernelModuleInlineTrans
59 
60  inline_trans = KernelModuleInlineTrans()
61  inline_trans.apply(schedule.walk(CodedKern)[0])
62 
63  print(schedule.parent.view())
64 
65 
66  .. warning ::
67  Not all kernel subroutines can be module-inlined. This transformation
68  will reject attempts to in-line kernels that access global data in the
69  original module.
70 
71  '''
72 
73  def __str__(self):
74  return "Inline a kernel subroutine into the PSy module"
75 
76  # pylint: disable=too-many-branches
77  def validate(self, node, options=None):
78  '''
79  Checks that the supplied node is a Kernel and that it is possible to
80  inline its PSyIR.
81 
82  :param kern: the kernel which is the target of the transformation.
83  :type kern: :py:class:`psyclone.psyGen.CodedKern`
84  :param options: a dictionary with options for transformations.
85  :type options: Optional[Dict[str, Any]]
86 
87  :raises TransformationError: if the target node is not a sub-class of \
88  psyGen.CodedKern.
89  :raises TransformationError: if the subroutine containing the \
90  implementation of the kernel cannot be retrieved with \
91  'get_kernel_schedule'.
92  :raises TransformationError: if the name of the routine that \
93  implements the kernel is not the same as the kernel name. This \
94  will happen if the kernel is polymorphic (uses a Fortran \
95  INTERFACE) and will be resolved by #1824.
96  :raises TransformationError: if the kernel cannot be safely inlined.
97 
98  '''
99  if not isinstance(node, CodedKern):
100  raise TransformationError(
101  f"Target of a {self.name} must be a sub-class of "
102  f"psyGen.CodedKern but got '{type(node).__name__}'")
103 
104  # Check that the PSyIR and associated Symbol table of the Kernel is OK.
105  # If this kernel contains symbols that are not captured in the PSyIR
106  # SymbolTable then this raises an exception.
107  try:
108  kernel_schedule = node.get_kernel_schedule()
109  except Exception as error:
110  raise TransformationError(
111  f"{self.name} failed to retrieve PSyIR for kernel "
112  f"'{node.name}' using the 'get_kernel_schedule' method"
113  f" due to {error}."
114  ) from error
115 
116  # We do not support kernels that use symbols representing global
117  # variables declared in its own parent module (we would need to
118  # create new imports to this module for those, and we don't do
119  # this yet).
120  # These can only be found in References, Calls and CodeBlocks
121  for var in kernel_schedule.walk(Reference):
122  symbol = var.symbol
123  if isinstance(symbol, IntrinsicSymbol):
124  continue
125  if not symbol.is_import:
126  try:
127  var.scope.symbol_table.lookup(
128  symbol.name, scope_limit=kernel_schedule)
129  except KeyError as err:
130  raise TransformationError(
131  f"Kernel '{node.name}' contains accesses to "
132  f"'{symbol.name}' which is declared in the same "
133  f"module scope. Cannot inline such a kernel.") from err
134  for block in kernel_schedule.walk(CodeBlock):
135  for name in block.get_symbol_names():
136  try:
137  block.scope.symbol_table.lookup(
138  name, scope_limit=kernel_schedule)
139  except KeyError as err:
140  if not block.scope.symbol_table.lookup(name).is_import:
141  raise TransformationError(
142  f"Kernel '{node.name}' contains accesses to "
143  f"'{name}' in a CodeBlock that is declared in the "
144  f"same module scope. "
145  f"Cannot inline such a kernel.") from err
146 
147  # We can't transform subroutines that shadow top-level symbol module
148  # names, because we won't be able to bring this into the subroutine
149  symtab = kernel_schedule.ancestor(Container).symbol_table
150  for scope in kernel_schedule.walk(ScopingNode):
151  for symbol in scope.symbol_table.symbols:
152  for mod in symtab.containersymbols:
153  if symbol.name == mod.name and not \
154  isinstance(symbol, ContainerSymbol):
155  raise TransformationError(
156  f"Kernel '{node.name}' cannot be module-inlined"
157  f" because the subroutine shadows the symbol "
158  f"name of the module container '{symbol.name}'.")
159 
160  # If the symbol already exist at the call site it must be referring
161  # to a Routine
162  try:
163  existing_symbol = node.scope.symbol_table.lookup(node.name)
164  except KeyError:
165  existing_symbol = None
166  if existing_symbol and not isinstance(existing_symbol, RoutineSymbol):
167  raise TransformationError(
168  f"Cannot module-inline subroutine '{node.name}' because "
169  f"symbol '{existing_symbol}' with the same name already "
170  f"exists and changing the name of module-inlined "
171  f"subroutines is not supported yet.")
172 
173  @staticmethod
174  def _prepare_code_to_inline(code_to_inline):
175  '''Prepare the PSyIR tree to inline by bringing in to the subroutine
176  all referenced symbols so that the implementation is self contained.
177 
178  TODO #2271 will improve this method and could potentially
179  avoid the need for debug_string() within get_kernel_schedule()
180  in dynamo0.3.py. Sergi suggests that we may be missing the
181  traversal of the declaration init expressions here and that
182  might solve the problem. I'm not so sure and explain why in
183  get_kernel_schedule() but still referencing this issue.
184 
185  :param code_to_inline: the subroutine to module-inline.
186  :type code_to_inline: :py:class:`psyclone.psyir.node.Routine`
187 
188  '''
189  # pylint: disable=too-many-branches
190  source_container = code_to_inline.ancestor(Container)
191 
192  # First make a set with all symbols used inside the subroutine
193  all_symbols = set()
194  for scope in code_to_inline.walk(ScopingNode):
195  for symbol in scope.symbol_table.symbols:
196  all_symbols.add(symbol)
197  for reference in code_to_inline.walk(Reference):
198  all_symbols.add(reference.symbol)
199  for literal in code_to_inline.walk(Literal):
200  # Literals may reference symbols in their precision
201  if isinstance(literal.datatype.precision, Symbol):
202  all_symbols.add(literal.datatype.precision)
203  for caller in code_to_inline.walk(Call):
204  all_symbols.add(caller.routine.symbol)
205  for cblock in code_to_inline.walk(CodeBlock):
206  for name in cblock.get_symbol_names():
207  all_symbols.add(cblock.scope.symbol_table.lookup(name))
208 
209  # Then decide which symbols need to be brought inside the subroutine
210  symbols_to_bring_in = set()
211  for symbol in all_symbols:
212  if symbol.is_unresolved or symbol.is_import:
213  # This symbol is already in the symbol table, but adding it
214  # to the 'symbols_to_bring_in' will make the next step bring
215  # into the subroutine all modules that it could come from.
216  symbols_to_bring_in.add(symbol)
217  if isinstance(symbol, DataSymbol):
218  # DataTypes can reference other symbols
219  if isinstance(symbol.datatype, DataTypeSymbol):
220  symbols_to_bring_in.add(symbol.datatype)
221  elif hasattr(symbol.datatype, 'precision'):
222  if isinstance(symbol.datatype.precision, Symbol):
223  symbols_to_bring_in.add(symbol.datatype.precision)
224 
225  # Bring the selected symbols inside the subroutine
226  for symbol in symbols_to_bring_in:
227  if symbol.name not in code_to_inline.symbol_table:
228  code_to_inline.symbol_table.add(symbol)
229  # And when necessary the modules where they come from
230  if symbol.is_unresolved:
231  # We don't know where this comes from, we need to bring
232  # in all top-level imports with wildcard imports
233  for mod in source_container.symbol_table.containersymbols:
234  if mod.wildcard_import:
235  if mod.name not in code_to_inline.symbol_table:
236  code_to_inline.symbol_table.add(mod)
237  else:
238  code_to_inline.symbol_table.lookup(mod.name).\
239  wildcard_import = True
240  elif symbol.is_import:
241  module_symbol = symbol.interface.container_symbol
242  if module_symbol.name not in code_to_inline.symbol_table:
243  code_to_inline.symbol_table.add(module_symbol)
244  else:
245  # If it already exists, we know its a container (from the
246  # validation) so we just need to point to it
247  symbol.interface.container_symbol = \
248  code_to_inline.symbol_table.lookup(module_symbol.name)
249 
250  def apply(self, node, options=None):
251  ''' Bring the kernel subroutine into this Container.
252 
253  :param node: the kernel to module-inline.
254  :type node: :py:class:`psyclone.psyGen.CodedKern`
255  :param options: a dictionary with options for transformations.
256  :type options: Optional[Dict[str, Any]]
257 
258  '''
259  self.validatevalidatevalidate(node, options)
260 
261  if not options:
262  options = {}
263 
264  # Note that we use the resolved callee subroutine name and not the
265  # caller one, this is important because if it is an interface it will
266  # use the concrete implementation name. When this happens the new name
267  # may already be in use, but the equality check below guarantees
268  # that if it exists it is only valid when it references the exact same
269  # implementation.
270  code_to_inline = node.get_kernel_schedule()
271  name = code_to_inline.name
272 
273  try:
274  existing_symbol = node.scope.symbol_table.lookup(name)
275  except KeyError:
276  existing_symbol = None
277 
278  self._prepare_code_to_inline_prepare_code_to_inline(code_to_inline)
279 
280  if not existing_symbol:
281  # If it doesn't exist already, module-inline the subroutine by:
282  # 1) Registering the subroutine symbol in the Container
283  node.ancestor(Container).symbol_table.add(RoutineSymbol(
284  name, interface=DefaultModuleInterface()
285  ))
286  # 2) Insert the relevant code into the tree.
287  node.ancestor(Container).addchild(code_to_inline.detach())
288  else:
289  # The routine symbol already exist, and we know from the validation
290  # that its a Routine. Now check if they are exactly the same.
291  for routine in node.ancestor(Container).walk(Routine,
292  stop_type=Routine):
293  if routine.name == node.name:
294  # This TransformationError happens here and not in the
295  # validation because it needs the symbols_to_bring_in
296  # applied to effectively compare both versions
297  # This will be fixed when module-inlining versioning is
298  # implemented.
299  if routine != code_to_inline:
300  raise TransformationError(
301  f"Cannot inline subroutine '{node.name}' because "
302  f"another, different, subroutine with the same "
303  f"name already exists and versioning of module-"
304  f"inlined subroutines is not implemented yet.")
305 
306  # We only modify the kernel call name after the equality check to
307  # ensure the apply will succeed and we don't leave with an inconsistent
308  # tree.
309  if node.name.lower() != name:
310  node.name = name
311 
312  # Set the module-inline flag to avoid generating the kernel imports
313  # TODO #1823. If the kernel imports were generated at PSy-layer
314  # creation time, we could just remove it here instead of setting a
315  # flag.
316  node.module_inline = True
def validate(self, node, options=None)
Definition: psyGen.py:2783