Reference Guide  2.5.0
lfric_kern.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 and S. Siso, STFC Daresbury Lab
35 # Modified I. Kavcic, A. Coughtrie, L. Turner and O. Brunt, Met Office
36 # Modified J. Henrichs, Bureau of Meteorology
37 # Modified A. B. G. Chalk and N. Nobre, STFC Daresbury Lab
38 
39 ''' This module implements the PSyclone LFRic API by specialising the required
40  base class Kern in psyGen.py '''
41 
42 # Imports
43 from collections import OrderedDict, namedtuple
44 
45 from psyclone.configuration import Config
46 from psyclone.core import AccessType
47 from psyclone.domain.lfric import (KernCallArgList, KernStubArgList,
48  KernelInterface, LFRicConstants)
49 from psyclone.errors import GenerationError, InternalError, FieldNotFoundError
50 from psyclone.f2pygen import ModuleGen, SubroutineGen, UseGen
51 from psyclone.parse.algorithm import Arg, KernelCall
52 from psyclone.psyGen import InvokeSchedule, CodedKern, args_filter
53 from psyclone.psyir.frontend.fparser2 import Fparser2Reader
54 from psyclone.psyir.nodes import (Loop, Literal, Reference,
55  KernelSchedule)
56 from psyclone.psyir.symbols import DataSymbol, ScalarType, ArrayType
57 
58 
60  ''' Stores information about LFRic Kernels as specified by the
61  Kernel metadata and associated algorithm call. Uses this
62  information to generate appropriate PSy layer code for the Kernel
63  instance or to generate a Kernel stub.
64 
65  '''
66  # pylint: disable=too-many-instance-attributes
67  # An instance of this `namedtuple` is used to store information on each of
68  # the quadrature rules required by a kernel.
69  #
70  # alg_name: The actual argument text specifying the QR object in the
71  # Alg. layer.
72  # psy_name: The PSy-layer variable name for the QR object.
73  # kernel_args: List of kernel arguments associated with this QR rule.
74 
75  QRRule = namedtuple("QRRule",
76  ["alg_name", "psy_name", "kernel_args"])
77 
78  def __init__(self):
79  # The super-init is called from the _setup() method which in turn
80  # is called from load().
81  # pylint: disable=super-init-not-called
82  # Import here to avoid circular dependency
83  # pylint: disable=import-outside-toplevel
84  if False: # pylint: disable=using-constant-test
85  from psyclone.dynamo0p3 import DynKernelArguments
86  self._arguments_arguments_arguments = DynKernelArguments(None, None) # for pyreverse
87  self._parent_parent = None
88  self._base_name_base_name = ""
89  self._func_descriptors_func_descriptors = None
90  self._fs_descriptors_fs_descriptors = None
91  # Whether this kernel requires quadrature
92  self._qr_required_qr_required = False
93  # Whether this kernel requires basis functions
94  self._basis_required_basis_required = False
95  # What shapes of evaluator/quadrature this kernel requires (if any)
96  self._eval_shapes_eval_shapes = []
97  # The function spaces on which to *evaluate* basis/diff-basis
98  # functions if an evaluator is required for this kernel. Is a dict with
99  # (mangled) FS names as keys and associated kernel argument as value.
100  self._eval_targets_eval_targets = OrderedDict()
101  # Will hold a dict of QRRule namedtuple objects, one for each QR
102  # rule required by a kernel, indexed by shape. Needs to be ordered
103  # because we must preserve the ordering specified in the metadata.
104  self._qr_rules_qr_rules = OrderedDict()
105  self._cma_operation_cma_operation = None
106  # Reference to the DynInterGrid object holding any inter-grid aspects
107  # of this kernel or None if it is not an intergrid kernel
108  self._intergrid_ref_intergrid_ref = None # Reference to this kernel inter-grid
109  # The reference-element properties required by this kernel
110  self._reference_element_reference_element = None
111  # The mesh properties required by this kernel
112  self._mesh_properties_mesh_properties = None
113  # Initialise kinds (precisions) of all kernel arguments (start
114  # with 'real' and 'integer' kinds)
115  api_config = Config.get().api_conf("dynamo0.3")
116  self._argument_kinds_argument_kinds = {api_config.default_kind["real"],
117  api_config.default_kind["integer"]}
118 
119  def reference_accesses(self, var_accesses):
120  '''Get all variable access information. All accesses are marked
121  according to the kernel metadata
122 
123  :param var_accesses: VariablesAccessInfo instance that stores the \
124  information about variable accesses.
125  :type var_accesses: \
126  :py:class:`psyclone.core.VariablesAccessInfo`
127 
128  '''
129  # Use the KernelCallArgList class, which can also provide variable
130  # access information:
131  create_arg_list = KernCallArgList(self)
132  create_arg_list.generate(var_accesses)
133 
134  super().reference_accesses(var_accesses)
135  # Set the current location index to the next location, since after
136  # this kernel a new statement starts.
137  var_accesses.next_location()
138 
139  def load(self, call, parent=None):
140  '''
141  Sets up kernel information with the call object which is
142  created by the parser. This object includes information about
143  the invoke call and the associated kernel.
144 
145  :param call: The KernelCall object from which to extract information
146  about this kernel
147  :type call: :py:class:`psyclone.parse.algorithm.KernelCall`
148  :param parent: The parent node of the kernel call in the AST
149  we are constructing. This will be a loop.
150  :type parent: :py:class:`psyclone.domain.lfric.LFRicLoop`
151  '''
152  self._setup_basis_setup_basis(call.ktype)
153  self._setup_setup(call.ktype, call.module_name, call.args, parent)
154 
155  def load_meta(self, ktype):
156  '''
157  Sets up kernel information with the kernel type object
158  which is created by the parser. The object includes the
159  metadata describing the kernel code.
160 
161  :param ktype: the kernel metadata object produced by the parser
162  :type ktype: :py:class:`psyclone.domain.lfric.LFRicKernMetadata`
163 
164  :raises InternalError: for an invalid data type of a scalar argument.
165  :raises GenerationError: if an invalid argument type is found \
166  in the kernel.
167 
168  '''
169  # pylint: disable=too-many-branches
170  # Create a name for each argument
171  args = []
172  const = LFRicConstants()
173  for idx, descriptor in enumerate(ktype.arg_descriptors):
174  pre = None
175  if descriptor.argument_type.lower() == "gh_operator":
176  pre = "op_"
177  elif descriptor.argument_type.lower() == "gh_columnwise_operator":
178  pre = "cma_op_"
179  elif descriptor.argument_type.lower() == "gh_field":
180  pre = "field_"
181  elif (descriptor.argument_type.lower() in
182  const.VALID_SCALAR_NAMES):
183  if descriptor.data_type.lower() == "gh_real":
184  pre = "rscalar_"
185  elif descriptor.data_type.lower() == "gh_integer":
186  pre = "iscalar_"
187  elif descriptor.data_type.lower() == "gh_logical":
188  pre = "lscalar_"
189  else:
190  raise InternalError(
191  f"Expected one of {const.VALID_SCALAR_DATA_TYPES} "
192  f"data types for a scalar argument but found "
193  f"'{descriptor.data_type}'.")
194  else:
195  raise GenerationError(
196  f"LFRicKern.load_meta() expected one of "
197  f"{const.VALID_ARG_TYPE_NAMES} but found "
198  f"'{descriptor.argument_type}'")
199  args.append(Arg("variable", pre+str(idx+1)))
200 
201  if descriptor.stencil:
202  if not descriptor.stencil["extent"]:
203  # Stencil size (in cells) is passed in
204  args.append(Arg("variable",
205  pre+str(idx+1)+"_stencil_size"))
206  if descriptor.stencil["type"] == "xory1d":
207  # Direction is passed in
208  args.append(Arg("variable", pre+str(idx+1)+"_direction"))
209 
210  # Initialise basis/diff basis so we can test whether quadrature
211  # or an evaluator is required
212  self._setup_basis_setup_basis(ktype)
213  if self._basis_required_basis_required:
214  for shape in self._eval_shapes_eval_shapes:
215  if shape in const.VALID_QUADRATURE_SHAPES:
216  # Add a quadrature argument for each required quadrature
217  # rule.
218  args.append(Arg("variable", "qr_"+shape))
219  self._setup_setup(ktype, "dummy_name", args, None, check=False)
220 
221  def _setup_basis(self, kmetadata):
222  '''
223  Initialisation of the basis/diff basis information. This may be
224  needed before general setup so is computed in a separate method.
225 
226  :param kmetadata: The kernel metadata object produced by the parser.
227  :type kmetadata: :py:class:`psyclone.domain.lfric.LFRicKernMetadata`
228  '''
229  for descriptor in kmetadata.func_descriptors:
230  if len(descriptor.operator_names) > 0:
231  self._basis_required_basis_required = True
232  self._eval_shapes_eval_shapes = kmetadata.eval_shapes[:]
233  break
234 
235  def _setup(self, ktype, module_name, args, parent, check=True):
236  # pylint: disable=too-many-arguments
237  # pylint: disable=too-many-branches, too-many-locals
238  '''Internal setup of kernel information.
239 
240  :param ktype: object holding information on the parsed metadata for \
241  this kernel.
242  :type ktype: :py:class:`psyclone.domain.lfric.LFRicKernMetadata`
243  :param str module_name: the name of the Fortran module that contains \
244  the source of this Kernel.
245  :param args: list of Arg objects produced by the parser for the \
246  arguments of this kernel call.
247  :type args: List[:py:class:`psyclone.parse.algorithm.Arg`]
248  :param parent: the parent of this kernel call in the generated \
249  AST (will be a loop object).
250  :type parent: :py:class:`psyclone.domain.lfric.LFRicLoop`
251  :param bool check: whether to check for consistency between the \
252  kernel metadata and the algorithm layer. Defaults to True.
253 
254  '''
255  # Import here to avoid circular dependency
256  # pylint: disable=import-outside-toplevel
257  from psyclone.dynamo0p3 import DynKernelArguments, FSDescriptors
258  super().__init__(DynKernelArguments,
259  KernelCall(module_name, ktype, args),
260  parent, check)
261  # Remove "_code" from the name if it exists to determine the
262  # base name which (if LFRic naming conventions are
263  # followed) is used as the root for the module and subroutine
264  # names.
265  if self.namenamenamename.lower().endswith("_code"):
266  self._base_name_base_name = self.namenamenamename[:-5]
267  else:
268  # TODO: #11 add a warning here when logging is added
269  self._base_name_base_name = self.namenamenamename
270  self._func_descriptors_func_descriptors = ktype.func_descriptors
271  # Keep a record of the type of CMA kernel identified when
272  # parsing the kernel metadata
273  self._cma_operation_cma_operation = ktype.cma_operation
274  self._fs_descriptors_fs_descriptors = FSDescriptors(ktype.func_descriptors)
275 
276  # If the kernel metadata specifies that this is an inter-grid kernel
277  # create the associated DynInterGrid
278  if ktype.is_intergrid:
279  if not self.ancestor(InvokeSchedule):
280  raise NotImplementedError(
281  f"Intergrid kernels can only be setup inside an "
282  f"InvokeSchedule, but attempted '{self.name}' without it.")
283  fine_args = args_filter(self.argumentsarguments.args,
284  arg_meshes=["gh_fine"])
285  coarse_args = args_filter(self.argumentsarguments.args,
286  arg_meshes=["gh_coarse"])
287 
288  from psyclone.dynamo0p3 import DynInterGrid
289  intergrid = DynInterGrid(fine_args[0], coarse_args[0])
290  self._intergrid_ref_intergrid_ref = intergrid
291 
292  const = LFRicConstants()
293  # Check that all specified evaluator shapes are recognised
294  invalid_shapes = set(self._eval_shapes_eval_shapes) \
295  - set(const.VALID_EVALUATOR_SHAPES)
296  if invalid_shapes:
297  raise InternalError(
298  f"Evaluator shape(s) {list(invalid_shapes)} is/are not "
299  f"recognised. Must be one of {const.VALID_EVALUATOR_SHAPES}.")
300 
301  # If there are any quadrature rule(s), what are the names of the
302  # corresponding algorithm arguments? Can't use set() here because
303  # we need to preserve the ordering specified in the metadata.
304  qr_shapes = [shape for shape in self._eval_shapes_eval_shapes if
305  shape in const.VALID_QUADRATURE_SHAPES]
306 
307  # The quadrature-related arguments to a kernel always come last so
308  # construct an enumerator with start value -<no. of qr rules>
309  for idx, shape in enumerate(qr_shapes, -len(qr_shapes)):
310 
311  qr_arg = args[idx]
312 
313  # Use the InvokeSchedule symbol_table to create a unique symbol
314  # name for the whole Invoke.
315  if qr_arg.varname:
316  tag = "AlgArgs_" + qr_arg.text
317  qr_name = self.ancestor(InvokeSchedule).symbol_table.\
318  find_or_create_integer_symbol(qr_arg.varname, tag=tag).name
319  else:
320  # If we don't have a name then we must be doing kernel-stub
321  # generation so create a suitable name.
322  # TODO #719 we don't yet have a symbol table to prevent
323  # clashes.
324  qr_name = "qr_"+shape.split("_")[-1]
325 
326  # LFRic api kernels require quadrature rule arguments to be
327  # passed in if one or more basis functions are used by the kernel
328  # and gh_shape == "gh_quadrature_***".
329  # if self._eval_shape == "gh_quadrature_xyz":
330  # self._qr_args = ["np_xyz", "weights_xyz"]
331  if shape == "gh_quadrature_xyoz":
332  qr_args = ["np_xy", "np_z", "weights_xy", "weights_z"]
333  # elif self._eval_shape == "gh_quadrature_xoyoz":
334  # qr_args = ["np_x", "np_y", "np_z",
335  # "weights_x", "weights_y", "weights_z"]
336  elif shape == "gh_quadrature_face":
337  qr_args = ["nfaces", "np_xyz", "weights_xyz"]
338  elif shape == "gh_quadrature_edge":
339  qr_args = ["nedges", "np_xyz", "weights_xyz"]
340  else:
341  raise InternalError(f"Unsupported quadrature shape "
342  f"('{shape}') found in LFRicKern._setup")
343 
344  # Append the name of the qr argument to the names of the qr-related
345  # variables.
346  qr_args = [arg + "_" + qr_name for arg in qr_args]
347 
348  self._qr_rules_qr_rules[shape] = self.QRRuleQRRule(qr_arg.text, qr_name, qr_args)
349 
350  if "gh_evaluator" in self._eval_shapes_eval_shapes:
351  # Kernel has an evaluator. If gh_evaluator_targets is present
352  # then that specifies the function spaces for which the evaluator
353  # is required. Otherwise, the FS of the updated argument(s) tells
354  # us upon which nodal points the evaluator will be required
355  for fs_name in ktype.eval_targets:
356  arg, fspace = self.argumentsarguments.get_arg_on_space_name(fs_name)
357  # Set up our dict of evaluator targets, one entry per
358  # target FS.
359  if fspace.mangled_name not in self._eval_targets_eval_targets:
360  self._eval_targets_eval_targets[fspace.mangled_name] = (fspace, arg)
361 
362  # Properties of the reference element required by this kernel
363  self._reference_element_reference_element = ktype.reference_element
364 
365  # Properties of the mesh required by this kernel
366  self._mesh_properties_mesh_properties = ktype.mesh
367 
368  @property
369  def qr_rules(self):
370  '''
371  :return: details of each of the quadrature rules required by this \
372  kernel.
373  :rtype: OrderedDict containing \
374  :py:class:`psyclone.domain.lfric.LFRicKern.QRRule` indexed by \
375  quadrature shape.
376  '''
377  return self._qr_rules_qr_rules
378 
379  @property
380  def cma_operation(self):
381  '''
382  :return: the type of CMA operation performed by this kernel
383  (one of 'assembly', 'apply' or 'matrix-matrix') or None
384  if the kernel does not involve CMA operators.
385  :rtype: str
386  '''
387  return self._cma_operation_cma_operation
388 
389  @property
390  def is_intergrid(self):
391  '''
392  :return: True if it is an inter-grid kernel, False otherwise
393  :rtype: bool
394  '''
395  return self._intergrid_ref_intergrid_ref is not None
396 
397  @property
398  def colourmap(self):
399  '''
400  Getter for the name of the colourmap associated with this kernel call.
401 
402  :returns: name of the colourmap (Fortran array).
403  :rtype: str
404 
405  :raises InternalError: if this kernel is not coloured.
406 
407  '''
408  if not self.is_colouredis_coloured():
409  raise InternalError(f"Kernel '{self.name}' is not inside a "
410  f"coloured loop.")
411  sched = self.ancestor(InvokeSchedule)
412  if self.is_intergridis_intergrid:
413  cmap = self._intergrid_ref_intergrid_ref.colourmap_symbol.name
414  else:
415  try:
416  cmap = sched.symbol_table.lookup_with_tag("cmap").name
417  except KeyError:
418  # We have to do this here as _init_colourmap (which calls this
419  # method) is only called at code-generation time.
420  cmap = sched.symbol_table.find_or_create_array(
421  "cmap", 2, ScalarType.Intrinsic.INTEGER,
422  tag="cmap").name
423 
424  return cmap
425 
426  @property
428  '''
429  Getter for the symbol of the array holding the index of the last
430  cell of each colour.
431 
432  :returns: name of the array.
433  :rtype: str
434 
435  :raises InternalError: if this kernel is not coloured.
436  '''
437  if not self.is_colouredis_coloured():
438  raise InternalError(f"Kernel '{self.name}' is not inside a "
439  f"coloured loop.")
440 
441  if self.is_intergridis_intergrid:
442  return self._intergrid_ref_intergrid_ref.last_cell_var_symbol
443 
444  ubnd_name = self.ancestor(Loop).upper_bound_name
445  const = LFRicConstants()
446 
447  if ubnd_name in const.HALO_ACCESS_LOOP_BOUNDS:
448  return self.scope.symbol_table.find_or_create_array(
449  "last_halo_cell_all_colours", 2,
450  ScalarType.Intrinsic.INTEGER,
451  tag="last_halo_cell_all_colours")
452 
453  return self.scope.symbol_table.find_or_create_array(
454  "last_edge_cell_all_colours", 1,
455  ScalarType.Intrinsic.INTEGER,
456  tag="last_edge_cell_all_colours")
457 
458  @property
459  def ncolours_var(self):
460  '''
461  Getter for the name of the variable holding the number of colours
462  associated with this kernel call.
463 
464  :return: name of the variable holding the number of colours
465  :rtype: Optional[str]
466 
467  :raises InternalError: if this kernel is not coloured.
468  '''
469  if not self.is_colouredis_coloured():
470  raise InternalError(f"Kernel '{self.name}' is not inside a "
471  f"coloured loop.")
472  if self.is_intergridis_intergrid:
473  ncols_sym = self._intergrid_ref_intergrid_ref.ncolours_var_symbol
474  if not ncols_sym:
475  return None
476  return ncols_sym.name
477 
478  return self.scope.symbol_table.lookup_with_tag("ncolour").name
479 
480  @property
481  def fs_descriptors(self):
482  '''
483  :return: a list of function space descriptor objects of
484  type FSDescriptors which contain information about
485  the function spaces.
486  :rtype: List[:py:class:`psyclone.FSDescriptors`].
487 
488  '''
489  return self._fs_descriptors_fs_descriptors
490 
491  @property
492  def qr_required(self):
493  '''
494  :return: True if this kernel requires quadrature, else returns False.
495  :rtype: bool
496 
497  '''
498  return self._basis_required_basis_required and self.qr_rulesqr_rules
499 
500  @property
501  def eval_shapes(self):
502  '''
503  :return: the value(s) of GH_SHAPE for this kernel or an empty list \
504  if none are specified.
505  :rtype: list
506 
507  '''
508  return self._eval_shapes_eval_shapes
509 
510  @property
511  def eval_targets(self):
512  '''
513  :return: the function spaces upon which basis/diff-basis functions \
514  are to be evaluated for this kernel.
515  :rtype: dict of (:py:class:`psyclone.domain.lfric.FunctionSpace`, \
516  :py:class`psyclone.dynamo0p3.DynKernelArgument`), indexed by \
517  the names of the target function spaces.
518  '''
519  return self._eval_targets_eval_targets
520 
521  @property
522  def reference_element(self):
523  '''
524  :returns: the reference-element properties required by this kernel.
525  :rtype: :py:class:`psyclone.dynamo0p3.RefElementMetaData`
526  '''
527  return self._reference_element_reference_element
528 
529  @property
530  def mesh(self):
531  '''
532  :returns: the mesh properties required by this kernel.
533  :rtype: :py:class`psyclone.dynamo0p3.MeshPropertiesMetaData`
534  '''
535  return self._mesh_properties_mesh_properties
536 
537  @property
539  '''
540  :returns: True if all of the arguments updated by this kernel have \
541  'GH_WRITE' access, False otherwise.
542  :rtype: bool
543 
544  '''
545  accesses = set(arg.access for arg in self.argsargs)
546  all_writes = AccessType.all_write_accesses()
547  all_writes.remove(AccessType.WRITE)
548  return (not accesses.intersection(set(all_writes)))
549 
550  @property
551  def base_name(self):
552  '''
553  :returns: a base name for this kernel.
554  :rtype: str
555  '''
556  return self._base_name_base_name
557 
558  @property
559  def argument_kinds(self):
560  '''
561  :returns: kinds (precisions) for all arguments in a kernel.
562  :rtype: set of str
563 
564  '''
565  return self._argument_kinds_argument_kinds
566 
567  @property
568  def gen_stub(self):
569  '''
570  Create the fparser1 AST for a kernel stub.
571 
572  :returns: root of fparser1 AST for the stub routine.
573  :rtype: :py:class:`fparser.one.block_statements.Module`
574 
575  :raises GenerationError: if the supplied kernel stub does not operate \
576  on a supported subset of the domain (currently only "cell_column").
577 
578  '''
579  # The operates-on/iterates-over values supported by the stub generator.
580  const = LFRicConstants()
581  supported_operates_on = const.USER_KERNEL_ITERATION_SPACES[:]
582  # TODO #925 Add support for 'domain' kernels
583  supported_operates_on.remove("domain")
584 
585  # Check operates-on (iteration space) before generating code
586  if self.iterates_overiterates_over not in supported_operates_on:
587  raise GenerationError(
588  f"The LFRic API kernel-stub generator supports kernels that "
589  f"operate on one of {supported_operates_on} but found "
590  f"'{self.iterates_over}' in kernel '{self.name}'.")
591 
592  # Create an empty PSy layer module
593  psy_module = ModuleGen(self._base_name_base_name+"_mod")
594 
595  # Create the subroutine
596  sub_stub = SubroutineGen(psy_module, name=self._base_name_base_name+"_code",
597  implicitnone=True)
598 
599  # Add all the declarations
600  # Import here to avoid circular dependency
601  # pylint: disable=import-outside-toplevel
602  from psyclone.domain.lfric import (LFRicScalarArgs, LFRicFields,
603  LFRicDofmaps, LFRicStencils)
604  from psyclone.dynamo0p3 import (DynCellIterators, DynFunctionSpaces,
605  DynCMAOperators, DynBoundaryConditions,
606  DynLMAOperators, LFRicMeshProperties,
607  DynBasisFunctions, DynReferenceElement)
608  for entities in [DynCellIterators, LFRicDofmaps, DynFunctionSpaces,
609  DynCMAOperators, LFRicScalarArgs, LFRicFields,
610  DynLMAOperators, LFRicStencils, DynBasisFunctions,
611  DynBoundaryConditions, DynReferenceElement,
612  LFRicMeshProperties]:
613  entities(self).declarations(sub_stub)
614 
615  # Add wildcard "use" statement for all supported argument
616  # kinds (precisions)
617  sub_stub.add(
618  UseGen(sub_stub,
619  name=const.UTILITIES_MOD_MAP["constants"]["module"]))
620 
621  # Create the arglist
622  create_arg_list = KernStubArgList(self)
623  create_arg_list.generate()
624 
625  # Add the arglist
626  sub_stub.args = create_arg_list.arglist
627 
628  # Add the subroutine to the parent module
629  psy_module.add(sub_stub)
630  return psy_module.root
631 
633  '''Returns a PSyIR Schedule representing the kernel code. The base
634  class creates the PSyIR schedule on first invocation which is
635  then checked for consistency with the kernel metadata
636  here. The Schedule is just generated on first invocation, this
637  allows us to retain transformations that may subsequently be
638  applied to the Schedule.
639 
640  Once issue #935 is implemented, this routine will return the
641  PSyIR Schedule using LFRic-specific PSyIR where possible.
642 
643  :returns: Schedule representing the kernel code.
644  :rtype: :py:class:`psyclone.psyGen.KernelSchedule`
645 
646  :raises GenerationError: if no subroutine matching this kernel can \
647  be found in the parse tree of the associated source code.
648  '''
649  if self._kern_schedule_kern_schedule_kern_schedule:
650  return self._kern_schedule_kern_schedule_kern_schedule
651 
652  # Get the PSyIR Kernel Schedule(s)
653  routines = Fparser2Reader().get_routine_schedules(self.namenamenamename, self.astast)
654  for routine in routines:
655  # If one of the symbols is not declared in a routine then
656  # this is only picked up when writing out the routine
657  # (raising a VisitorError), so we check here so that
658  # invalid code is not inlined. We use debug_string() to
659  # minimise the overhead.
660 
661  # TODO #2271 could potentially avoid the need for
662  # debug_string() within. Sergi suggests that we may be
663  # missing the traversal of the declaration init
664  # expressions and that might solve the problem. I'm not so
665  # sure as we are talking about unknown symbols that will
666  # only be resolved in the back-end (or not). If I am right
667  # then one option would be to use the FortranWriter, but
668  # that would be bigger overhead, or perhaps just the
669  # declarations part of FortranWriter if that is possible.
670  # Also see TODO issue #2336 which captures the specific
671  # problem in LFRic that this fixes.
672  routine.debug_string()
673 
674  if len(routines) == 1:
675  sched = routines[0]
676  # TODO #928: We don't validate the arguments yet because the
677  # validation has many false negatives.
678  # self.validate_kernel_code_args(sched.symbol_table)
679  else:
680  # The kernel name corresponds to an interface block. Find which
681  # of the routines matches the precision of the arguments.
682  for routine in routines:
683  try:
684  # The validity check for the kernel arguments will raise
685  # an exception if the precisions don't match.
686  self.validate_kernel_code_argsvalidate_kernel_code_args(routine.symbol_table)
687  sched = routine
688  break
689  except GenerationError:
690  pass
691  else:
692  raise GenerationError(
693  f"Failed to find a kernel implementation with an interface"
694  f" that matches the invoke of '{self.name}'. (Tried "
695  f"routines {[item.name for item in routines]}.)")
696 
697  # TODO #935 - replace the PSyIR argument data symbols with LFRic data
698  # symbols. For the moment we just return the unmodified PSyIR schedule
699  # but this should use RaisePSyIR2LFRicKernTrans once KernelInterface
700  # is fully functional (#928).
701  ksched = KernelSchedule(sched.name,
702  symbol_table=sched.symbol_table.detach())
703  for child in sched.pop_all_children():
704  ksched.addchild(child)
705  sched.replace_with(ksched)
706 
707  self._kern_schedule_kern_schedule_kern_schedule = ksched
708 
709  return self._kern_schedule_kern_schedule_kern_schedule
710 
711  def validate_kernel_code_args(self, table):
712  '''Check that the arguments in the kernel code match the expected
713  arguments as defined by the kernel metadata and the LFRic
714  API.
715 
716  :param table: the symbol table to validate against the metadata.
717  :type table: :py:class:`psyclone.psyir.symbols.SymbolTable`
718 
719  :raises GenerationError: if the number of arguments indicated by the \
720  kernel metadata doesn't match the actual number of arguments in \
721  the symbol table.
722 
723  '''
724  # Get the kernel subroutine arguments
725  kern_code_args = table.argument_list
726 
727  # Get the kernel code interface according to the kernel
728  # metadata and LFRic API
729  interface_info = KernelInterface(self)
730  interface_info.generate()
731  interface_args = interface_info.arglist
732 
733  # 1: Check the the number of arguments match
734  actual_n_args = len(kern_code_args)
735  expected_n_args = len(interface_args)
736  if actual_n_args != expected_n_args:
737  raise GenerationError(
738  f"In kernel '{self.name}' the number of arguments indicated by"
739  f" the kernel metadata is {expected_n_args} but the actual "
740  f"number of kernel arguments found is {actual_n_args}.")
741 
742  # 2: Check that the properties of each argument match.
743  for idx, kern_code_arg in enumerate(kern_code_args):
744  interface_arg = interface_args[idx]
745  try:
746  alg_idx = interface_info.metadata_index_from_actual_index(idx)
747  alg_arg = self.argumentsarguments.args[alg_idx]
748  except KeyError:
749  # There's no algorithm argument directly associated with this
750  # kernel argument. (We only care about the data associated
751  # with scalar, field and operator arguments.)
752  alg_arg = None
753  self._validate_kernel_code_arg_validate_kernel_code_arg(kern_code_arg, interface_arg,
754  alg_arg)
755 
756  def _validate_kernel_code_arg(self, kern_code_arg, interface_arg,
757  alg_arg=None):
758  '''Internal method to check that the supplied argument descriptions
759  match and raise appropriate exceptions if not.
760 
761  :param kern_code_arg: kernel code argument.
762  :type kern_code_arg: :py:class:`psyclone.psyir.symbols.DataSymbol`
763  :param interface_arg: expected argument.
764  :type interface_arg: :py:class:`psyclone.psyir.symbols.DataSymbol`
765  :param alg_arg: the associated argument in the Algorithm layer. Note \
766  that only kernel arguments holding the data associated with \
767  scalar, field and operator arguments directly correspond to \
768  arguments that appear in the Algorithm layer.
769  :type alg_arg: \
770  Optional[:py:class`psyclone.dynamo0p3.DynKernelArgument`]
771 
772  :raises GenerationError: if the contents of the arguments do \
773  not match.
774  :raises InternalError: if an unexpected datatype is found.
775 
776  '''
777  # 1: intrinsic datatype
778  actual_datatype = kern_code_arg.datatype.intrinsic
779  expected_datatype = interface_arg.datatype.intrinsic
780  if actual_datatype != expected_datatype:
781  raise GenerationError(
782  f"Kernel argument '{kern_code_arg.name}' has datatype "
783  f"'{actual_datatype}' in kernel '{self.name}' but the LFRic "
784  f"API expects '{expected_datatype}'.")
785  # 2: precision. An LFRic kernel is only permitted to have a precision
786  # specified by a recognised type parameter or a no. of bytes.
787  actual_precision = kern_code_arg.datatype.precision
788  api_config = Config.get().api_conf("dynamo0.3")
789  if isinstance(actual_precision, DataSymbol):
790  # Convert precision into number of bytes to support
791  # mixed-precision kernels.
792  # TODO #1941: it would be better if the LFRic constants_mod.f90
793  # was the single source of truth for precision values.
794  actual_precision = api_config.precision_map[
795  actual_precision.name]
796  elif not isinstance(actual_precision, int):
797  raise GenerationError(
798  f"An argument to an LFRic kernel must have a precision defined"
799  f" by either a recognised LFRic type parameter (one of "
800  f"{sorted(api_config.precision_map.keys())}) or an integer"
801  f" number of bytes but argument '{kern_code_arg.name}' to "
802  f"kernel '{self.name}' has precision {actual_precision}.")
803 
804  if alg_arg:
805  # We have information on the corresponding argument in the
806  # Algorithm layer so we can check that the precision matches.
807  # This is used to identify the correct kernel subroutine for a
808  # mixed-precision kernel.
809  alg_precision = api_config.precision_map[alg_arg.precision]
810  if alg_precision != actual_precision:
811  raise GenerationError(
812  f"Precision ({alg_precision} bytes) of algorithm-layer "
813  f"argument '{alg_arg.name}' does not match that "
814  f"({actual_precision} bytes) of the corresponding kernel "
815  f"subroutine argument '{kern_code_arg.name}' for kernel "
816  f"'{self.name}'.")
817 
818  # 3: intent
819  actual_intent = kern_code_arg.interface.access
820  expected_intent = interface_arg.interface.access
821  if actual_intent.name != expected_intent.name:
822  raise GenerationError(
823  f"Kernel argument '{kern_code_arg.name}' has intent "
824  f"'{actual_intent.name}' in kernel '{self.name}' but the "
825  f"LFRic API expects intent '{expected_intent.name}'.")
826  # 4: scalar or array
827  if interface_arg.is_scalar:
828  if not kern_code_arg.is_scalar:
829  raise GenerationError(
830  f"Argument '{kern_code_arg.name}' to kernel '{self.name}' "
831  f"should be a scalar according to the LFRic API, but it "
832  f"is not.")
833  elif interface_arg.is_array:
834  if not kern_code_arg.is_array:
835  raise GenerationError(
836  f"Argument '{kern_code_arg.name}' to kernel '{self.name}' "
837  f"should be an array according to the LFRic API, but it "
838  f"is not.")
839  # 4.1: array dimensions
840  if len(interface_arg.shape) != len(kern_code_arg.shape):
841  raise GenerationError(
842  f"Argument '{kern_code_arg.name}' to kernel '{self.name}' "
843  f"should be an array with {len(interface_arg.shape)} "
844  f"dimension(s) according to the LFRic API, but "
845  f"found {len(kern_code_arg.shape)}.")
846  for dim_idx, kern_code_arg_dim in enumerate(kern_code_arg.shape):
847  if not isinstance(kern_code_arg_dim, ArrayType.ArrayBounds):
848  continue
849  if (not isinstance(kern_code_arg_dim.lower, Literal) or
850  kern_code_arg_dim.lower.value != "1"):
851  raise GenerationError(
852  f"All array arguments to LFRic kernels must have lower"
853  f" bounds of 1 for all dimensions. However, array "
854  f"'{kern_code_arg.name}' has a lower bound of "
855  f"'{kern_code_arg_dim.lower}' for dimension {dim_idx}")
856  kern_code_arg_upper_dim = kern_code_arg_dim.upper
857  interface_arg_upper_dim = interface_arg.shape[dim_idx].upper
858  if (isinstance(kern_code_arg_upper_dim, Reference) and
859  isinstance(interface_arg_upper_dim, Reference) and
860  isinstance(kern_code_arg_upper_dim.symbol,
861  DataSymbol) and
862  isinstance(interface_arg_upper_dim.symbol,
863  DataSymbol)):
864  # Only check when there is a symbol. Unspecified
865  # dimensions, dimensions with scalar values,
866  # offsets, or dimensions that include arithmetic
867  # are skipped.
868  try:
869  self._validate_kernel_code_arg_validate_kernel_code_arg(
870  kern_code_arg_upper_dim.symbol,
871  interface_arg_upper_dim.symbol)
872  except GenerationError as info:
873  raise GenerationError(
874  f"For dimension {dim_idx+1} in array argument "
875  f"'{kern_code_arg.name}' to kernel '{self.name}' "
876  f"the following error was found: "
877  f"{info.args[0]}") from info
878  else:
879  raise InternalError(
880  f"Unexpected argument type found for '{kern_code_arg.name}' in"
881  f" kernel '{self.name}'. Expecting a scalar or an array.")
882 
884  '''
885  Perform validation checks for any global constraints (that require the
886  tree to be complete).
887 
888  :raises GenerationError: if this kernel does not have a supported
889  operates-on (currently only "cell_column").
890  :raises GenerationError: if the loop goes beyond the level 1
891  halo and an operator is accessed.
892  :raises GenerationError: if a kernel in the loop has an inc access
893  and the loop is not coloured but is within an OpenMP
894  parallel region.
895  '''
896  # Check operates-on (iteration space) before generating code
897  const = LFRicConstants()
898  if self.iterates_overiterates_over not in const.USER_KERNEL_ITERATION_SPACES:
899  raise GenerationError(
900  f"The LFRic API supports calls to user-supplied kernels that "
901  f"operate on one of {const.USER_KERNEL_ITERATION_SPACES}, but "
902  f"kernel '{self.name}' operates on '{self.iterates_over}'.")
903 
904  # pylint: disable=import-outside-toplevel
905  from psyclone.domain.lfric import LFRicLoop
906  parent_loop = self.ancestor(LFRicLoop)
907 
908  # Check whether this kernel reads from an operator
909  op_args = parent_loop.args_filter(
910  arg_types=const.VALID_OPERATOR_NAMES,
911  arg_accesses=[AccessType.READ, AccessType.READWRITE])
912  if op_args:
913  # It does. We must check that our parent loop does not
914  # go beyond the L1 halo.
915  if (parent_loop.upper_bound_name == "cell_halo" and
916  parent_loop.upper_bound_halo_depth > 1):
917  raise GenerationError(
918  f"Kernel '{self._name}' reads from an operator and "
919  f"therefore cannot be used for cells beyond the level 1 "
920  f"halo. However the containing loop goes out to level "
921  f"{parent_loop.upper_bound_halo_depth}")
922 
923  if not self.is_colouredis_coloured():
924  # This kernel call has not been coloured
925  # - is it OpenMP parallel, i.e. are we a child of
926  # an OpenMP directive?
927  if self.is_openmp_parallel():
928  try:
929  # It is OpenMP parallel - does it have an argument
930  # with INC access?
931  _ = self.incremented_argincremented_arg()
932  raise GenerationError(f"Kernel '{self._name}' has an "
933  f"argument with INC access and "
934  f"therefore must be coloured in "
935  f"order to be parallelised with "
936  f"OpenMP.")
937  except FieldNotFoundError:
938  pass
939 
941 
942 
943 # ---------- Documentation utils -------------------------------------------- #
944 # The list of module members that we wish AutoAPI to generate
945 # documentation for. (See https://psyclone-ref.readthedocs.io)
946 __all__ = ['LFRicKern']
def _validate_kernel_code_arg(self, kern_code_arg, interface_arg, alg_arg=None)
Definition: lfric_kern.py:757
def load(self, call, parent=None)
Definition: lfric_kern.py:139
def _setup(self, ktype, module_name, args, parent, check=True)
Definition: lfric_kern.py:235
def reference_accesses(self, var_accesses)
Definition: lfric_kern.py:119
def incremented_arg(self)
Definition: psyGen.py:1567
def name(self, value)
Definition: psyGen.py:1321
def name(self)
Definition: psyGen.py:1313
def arguments(self)
Definition: psyGen.py:1309
def args(self)
Definition: psyGen.py:1080
def is_coloured(self)
Definition: psyGen.py:1329
def iterates_over(self)
Definition: psyGen.py:1343