39 ''' This module implements the PSyclone LFRic API by specialising the required
40 base class Kern in psyGen.py '''
43 from collections
import OrderedDict, namedtuple
48 KernelInterface, LFRicConstants)
49 from psyclone.errors import GenerationError, InternalError, FieldNotFoundError
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.
75 QRRule = namedtuple(
"QRRule",
76 [
"alg_name",
"psy_name",
"kernel_args"])
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"]}
120 '''Get all variable access information. All accesses are marked
121 according to the kernel metadata
123 :param var_accesses: VariablesAccessInfo instance that stores the \
124 information about variable accesses.
125 :type var_accesses: \
126 :py:class:`psyclone.core.VariablesAccessInfo`
132 create_arg_list.generate(var_accesses)
137 var_accesses.next_location()
139 def load(self, call, parent=None):
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.
145 :param call: The KernelCall object from which to extract information
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`
153 self.
_setup_setup(call.ktype, call.module_name, call.args, parent)
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.
161 :param ktype: the kernel metadata object produced by the parser
162 :type ktype: :py:class:`psyclone.domain.lfric.LFRicKernMetadata`
164 :raises InternalError: for an invalid data type of a scalar argument.
165 :raises GenerationError: if an invalid argument type is found \
173 for idx, descriptor
in enumerate(ktype.arg_descriptors):
175 if descriptor.argument_type.lower() ==
"gh_operator":
177 elif descriptor.argument_type.lower() ==
"gh_columnwise_operator":
179 elif descriptor.argument_type.lower() ==
"gh_field":
181 elif (descriptor.argument_type.lower()
in
182 const.VALID_SCALAR_NAMES):
183 if descriptor.data_type.lower() ==
"gh_real":
185 elif descriptor.data_type.lower() ==
"gh_integer":
187 elif descriptor.data_type.lower() ==
"gh_logical":
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}'.")
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)))
201 if descriptor.stencil:
202 if not descriptor.stencil[
"extent"]:
204 args.append(
Arg(
"variable",
205 pre+str(idx+1)+
"_stencil_size"))
206 if descriptor.stencil[
"type"] ==
"xory1d":
208 args.append(
Arg(
"variable", pre+str(idx+1)+
"_direction"))
215 if shape
in const.VALID_QUADRATURE_SHAPES:
218 args.append(
Arg(
"variable",
"qr_"+shape))
219 self.
_setup_setup(ktype,
"dummy_name", args,
None, check=
False)
221 def _setup_basis(self, kmetadata):
223 Initialisation of the basis/diff basis information. This may be
224 needed before general setup so is computed in a separate method.
226 :param kmetadata: The kernel metadata object produced by the parser.
227 :type kmetadata: :py:class:`psyclone.domain.lfric.LFRicKernMetadata`
229 for descriptor
in kmetadata.func_descriptors:
230 if len(descriptor.operator_names) > 0:
232 self.
_eval_shapes_eval_shapes = kmetadata.eval_shapes[:]
235 def _setup(self, ktype, module_name, args, parent, check=True):
238 '''Internal setup of kernel information.
240 :param ktype: object holding information on the parsed metadata for \
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.
258 super().__init__(DynKernelArguments,
274 self.
_fs_descriptors_fs_descriptors = FSDescriptors(ktype.func_descriptors)
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"])
289 intergrid = DynInterGrid(fine_args[0], coarse_args[0])
295 - set(const.VALID_EVALUATOR_SHAPES)
298 f
"Evaluator shape(s) {list(invalid_shapes)} is/are not "
299 f
"recognised. Must be one of {const.VALID_EVALUATOR_SHAPES}.")
304 qr_shapes = [shape
for shape
in self.
_eval_shapes_eval_shapes
if
305 shape
in const.VALID_QUADRATURE_SHAPES]
309 for idx, shape
in enumerate(qr_shapes, -len(qr_shapes)):
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
324 qr_name =
"qr_"+shape.split(
"_")[-1]
331 if shape ==
"gh_quadrature_xyoz":
332 qr_args = [
"np_xy",
"np_z",
"weights_xy",
"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"]
342 f
"('{shape}') found in LFRicKern._setup")
346 qr_args = [arg +
"_" + qr_name
for arg
in qr_args]
348 self.
_qr_rules_qr_rules[shape] = self.
QRRuleQRRule(qr_arg.text, qr_name, qr_args)
355 for fs_name
in ktype.eval_targets:
356 arg, fspace = self.
argumentsarguments.get_arg_on_space_name(fs_name)
359 if fspace.mangled_name
not in self.
_eval_targets_eval_targets:
360 self.
_eval_targets_eval_targets[fspace.mangled_name] = (fspace, arg)
371 :return: details of each of the quadrature rules required by this \
373 :rtype: OrderedDict containing \
374 :py:class:`psyclone.domain.lfric.LFRicKern.QRRule` indexed by \
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.
392 :return: True if it is an inter-grid kernel, False otherwise
400 Getter for the name of the colourmap associated with this kernel call.
402 :returns: name of the colourmap (Fortran array).
405 :raises InternalError: if this kernel is not coloured.
411 sched = self.ancestor(InvokeSchedule)
416 cmap = sched.symbol_table.lookup_with_tag(
"cmap").name
420 cmap = sched.symbol_table.find_or_create_array(
421 "cmap", 2, ScalarType.Intrinsic.INTEGER,
429 Getter for the symbol of the array holding the index of the last
432 :returns: name of the array.
435 :raises InternalError: if this kernel is not coloured.
444 ubnd_name = self.ancestor(Loop).upper_bound_name
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")
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")
461 Getter for the name of the variable holding the number of colours
462 associated with this kernel call.
464 :return: name of the variable holding the number of colours
465 :rtype: Optional[str]
467 :raises InternalError: if this kernel is not coloured.
476 return ncols_sym.name
478 return self.scope.symbol_table.lookup_with_tag(
"ncolour").name
483 :return: a list of function space descriptor objects of
484 type FSDescriptors which contain information about
486 :rtype: List[:py:class:`psyclone.FSDescriptors`].
494 :return: True if this kernel requires quadrature, else returns False.
503 :return: the value(s) of GH_SHAPE for this kernel or an empty list \
504 if none are specified.
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.
524 :returns: the reference-element properties required by this kernel.
525 :rtype: :py:class:`psyclone.dynamo0p3.RefElementMetaData`
532 :returns: the mesh properties required by this kernel.
533 :rtype: :py:class`psyclone.dynamo0p3.MeshPropertiesMetaData`
540 :returns: True if all of the arguments updated by this kernel have \
541 'GH_WRITE' access, False otherwise.
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)))
553 :returns: a base name for this kernel.
561 :returns: kinds (precisions) for all arguments in a kernel.
570 Create the fparser1 AST for a kernel stub.
572 :returns: root of fparser1 AST for the stub routine.
573 :rtype: :py:class:`fparser.one.block_statements.Module`
575 :raises GenerationError: if the supplied kernel stub does not operate \
576 on a supported subset of the domain (currently only "cell_column").
581 supported_operates_on = const.USER_KERNEL_ITERATION_SPACES[:]
583 supported_operates_on.remove(
"domain")
586 if self.
iterates_overiterates_over
not in supported_operates_on:
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}'.")
603 LFRicDofmaps, LFRicStencils)
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)
619 name=const.UTILITIES_MOD_MAP[
"constants"][
"module"]))
623 create_arg_list.generate()
626 sub_stub.args = create_arg_list.arglist
629 psy_module.add(sub_stub)
630 return psy_module.root
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.
640 Once issue #935 is implemented, this routine will return the
641 PSyIR Schedule using LFRic-specific PSyIR where possible.
643 :returns: Schedule representing the kernel code.
644 :rtype: :py:class:`psyclone.psyGen.KernelSchedule`
646 :raises GenerationError: if no subroutine matching this kernel can \
647 be found in the parse tree of the associated source code.
654 for routine
in routines:
672 routine.debug_string()
674 if len(routines) == 1:
682 for routine
in routines:
689 except 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]}.)")
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)
712 '''Check that the arguments in the kernel code match the expected
713 arguments as defined by the kernel metadata and the LFRic
716 :param table: the symbol table to validate against the metadata.
717 :type table: :py:class:`psyclone.psyir.symbols.SymbolTable`
719 :raises GenerationError: if the number of arguments indicated by the \
720 kernel metadata doesn't match the actual number of arguments in \
725 kern_code_args = table.argument_list
730 interface_info.generate()
731 interface_args = interface_info.arglist
734 actual_n_args = len(kern_code_args)
735 expected_n_args = len(interface_args)
736 if actual_n_args != expected_n_args:
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}.")
743 for idx, kern_code_arg
in enumerate(kern_code_args):
744 interface_arg = interface_args[idx]
746 alg_idx = interface_info.metadata_index_from_actual_index(idx)
747 alg_arg = self.
argumentsarguments.args[alg_idx]
756 def _validate_kernel_code_arg(self, kern_code_arg, interface_arg,
758 '''Internal method to check that the supplied argument descriptions
759 match and raise appropriate exceptions if not.
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.
770 Optional[:py:class`psyclone.dynamo0p3.DynKernelArgument`]
772 :raises GenerationError: if the contents of the arguments do \
774 :raises InternalError: if an unexpected datatype is found.
778 actual_datatype = kern_code_arg.datatype.intrinsic
779 expected_datatype = interface_arg.datatype.intrinsic
780 if actual_datatype != expected_datatype:
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}'.")
787 actual_precision = kern_code_arg.datatype.precision
788 api_config = Config.get().api_conf(
"dynamo0.3")
789 if isinstance(actual_precision, DataSymbol):
794 actual_precision = api_config.precision_map[
795 actual_precision.name]
796 elif not isinstance(actual_precision, int):
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}.")
809 alg_precision = api_config.precision_map[alg_arg.precision]
810 if alg_precision != actual_precision:
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 "
819 actual_intent = kern_code_arg.interface.access
820 expected_intent = interface_arg.interface.access
821 if actual_intent.name != expected_intent.name:
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}'.")
827 if interface_arg.is_scalar:
828 if not kern_code_arg.is_scalar:
830 f
"Argument '{kern_code_arg.name}' to kernel '{self.name}' "
831 f
"should be a scalar according to the LFRic API, but it "
833 elif interface_arg.is_array:
834 if not kern_code_arg.is_array:
836 f
"Argument '{kern_code_arg.name}' to kernel '{self.name}' "
837 f
"should be an array according to the LFRic API, but it "
840 if len(interface_arg.shape) != len(kern_code_arg.shape):
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):
849 if (
not isinstance(kern_code_arg_dim.lower, Literal)
or
850 kern_code_arg_dim.lower.value !=
"1"):
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,
862 isinstance(interface_arg_upper_dim.symbol,
870 kern_code_arg_upper_dim.symbol,
871 interface_arg_upper_dim.symbol)
872 except GenerationError
as info:
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
880 f
"Unexpected argument type found for '{kern_code_arg.name}' in"
881 f
" kernel '{self.name}'. Expecting a scalar or an array.")
885 Perform validation checks for any global constraints (that require the
886 tree to be complete).
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
898 if self.
iterates_overiterates_over
not in const.USER_KERNEL_ITERATION_SPACES:
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}'.")
906 parent_loop = self.ancestor(LFRicLoop)
909 op_args = parent_loop.args_filter(
910 arg_types=const.VALID_OPERATOR_NAMES,
911 arg_accesses=[AccessType.READ, AccessType.READWRITE])
915 if (parent_loop.upper_bound_name ==
"cell_halo" and
916 parent_loop.upper_bound_halo_depth > 1):
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}")
927 if self.is_openmp_parallel():
933 f
"argument with INC access and "
934 f
"therefore must be coloured in "
935 f
"order to be parallelised with "
937 except FieldNotFoundError:
946 __all__ = [
'LFRicKern']
def _validate_kernel_code_arg(self, kern_code_arg, interface_arg, alg_arg=None)
def load(self, call, parent=None)
def get_kernel_schedule(self)
def load_meta(self, ktype)
def last_cell_all_colours_symbol(self)
def validate_global_constraints(self)
def _setup_basis(self, kmetadata)
def reference_element(self)
def validate_kernel_code_args(self, table)
def _setup(self, ktype, module_name, args, parent, check=True)
def reference_accesses(self, var_accesses)
def all_updates_are_writes(self)
def incremented_arg(self)