Reference Guide  2.5.0
kern_call_arg_list.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 and L. Turner, Met Office
36 # Modified J. Henrichs, Bureau of Meteorology
37 
38 '''This module implements a class that manages the argument for a kernel
39 call. It especially adds all implicitly required parameters.
40 It creates the argument in two formats: first as a list of strings, but also
41 as a list of PSyIR nodes. TODO #1930: the support for the string format
42 should be removed as we migrate to use PSyIR in LFRic.
43 '''
44 
45 from collections import namedtuple
46 
47 from psyclone import psyGen
48 from psyclone.core import AccessType, Signature
49 from psyclone.domain.lfric import ArgOrdering, LFRicConstants
50 # Avoid circular import:
51 from psyclone.domain.lfric.lfric_types import LFRicTypes
52 from psyclone.errors import GenerationError, InternalError
53 from psyclone.psyir.nodes import ArrayReference, Reference, StructureReference
54 from psyclone.psyir.symbols import (
55  DataSymbol, DataTypeSymbol, UnresolvedType, ContainerSymbol,
56  ImportInterface, ScalarType)
57 
58 # psyir has classes created at runtime
59 # pylint: disable=no-member
60 # pylint: disable=too-many-lines
61 
62 
64  # pylint: disable=too-many-public-methods
65  # TODO: #845 Check that all implicit variables have the right type.
66  '''Creates the argument list required to call kernel "kern" from the
67  PSy-layer and captures the positions of the following arguments in
68  the argument list: nlayers, number of quadrature points and number
69  of degrees of freedom. The ordering and type of arguments is
70  captured by the base class.
71 
72  :param kern: The kernel that is being called.
73  :type kern: :py:class:`psyclone.domain.lfric.LFRicKern`
74 
75  '''
76  NdfInfo = namedtuple("NdfInfo", ["position", "function_space"])
77 
78  def __init__(self, kern):
79  super().__init__(kern)
80  self._nlayers_positions_nlayers_positions = []
81  self._nqp_positions_nqp_positions = []
82  self._ndf_positions_ndf_positions = []
83 
84  def get_user_type(self, module_name, user_type, name, tag=None):
85  # pylint: disable=too-many-arguments
86  '''Returns the symbol for a user-defined type. If required, the
87  required import statements will all be generated.
88 
89  :param str module_name: the name of the module from which the \
90  user-defined type must be imported.
91  :param str user_type: the name of the user-defined type.
92  :param str name: the name of the variable to be used in the Reference.
93  :param Optional[str] tag: tag to use for the variable, defaults to \
94  the name
95 
96  :return: the symbol that is used in the reference
97  :rtype: :py:class:`psyclone.psyir.symbols.Symbol`
98 
99  '''
100  if not tag:
101  tag = name
102 
103  try:
104  sym = self._symtab_symtab.lookup_with_tag(tag)
105  return sym
106  except KeyError:
107  pass
108 
109  # The symbol does not exist already. So we potentially need to
110  # create the import statement for the type:
111  try:
112  # Check if the module is already declared:
113  module = self._symtab_symtab.lookup(module_name)
114  except KeyError:
115  module = self._symtab_symtab.new_symbol(module_name,
116  symbol_type=ContainerSymbol)
117 
118  # Get the symbol table in which the module is declared:
119  mod_sym_tab = module.find_symbol_table(self._kern_kern)
120 
121  # The user-defined type must be declared in the same symbol
122  # table as the container (otherwise errors will happen later):
123  user_type_symbol = mod_sym_tab.find_or_create(
124  user_type,
125  symbol_type=DataTypeSymbol,
126  datatype=UnresolvedType(),
127  interface=ImportInterface(module))
128  # Declare the actual user symbol in the local symbol table, using
129  # the datatype from the root table:
130  sym = self._symtab_symtab.new_symbol(name, tag=tag,
131  symbol_type=DataSymbol,
132  datatype=user_type_symbol)
133  return sym
134 
135  def append_structure_reference(self, module_name, user_type, member_list,
136  name, tag=None, overwrite_datatype=None):
137  # pylint: disable=too-many-arguments
138  '''Creates a reference to a variable of a user-defined type. If
139  required, the required import statements will all be generated.
140 
141  :param str module_name: the name of the module from which the
142  user-defined type must be imported.
143  :param str user_type: the name of the user-defined type.
144  :param member_list: the members used hierarchically.
145  :type member_list: List[str]
146  :param str name: the name of the variable to be used in the Reference.
147  :param Optional[str] tag: tag to use for the variable, defaults to
148  the name
149  :param overwrite_datatype: the datatype for the reference, which will
150  overwrite the value determined by analysing the corresponding
151  user defined type. This is useful when e.g. the module that
152  declares the structure cannot be accessed.
153  :type overwrite_datatype:
154  Optional[:py:class:`psyclone.psyir.symbols.DataType`]
155 
156  :return: the symbol that is used in the reference
157  :rtype: :py:class:`psyclone.psyir.symbols.Symbol`
158 
159  '''
160  sym = self.get_user_typeget_user_type(module_name, user_type, name, tag)
161  self.psyir_appendpsyir_append(StructureReference.
162  create(sym, member_list,
163  overwrite_datatype=overwrite_datatype))
164  return sym
165 
166  def cell_position(self, var_accesses=None):
167  '''Adds a cell argument to the argument list and if supplied stores
168  this access in var_accesses.
169 
170  :param var_accesses: optional VariablesAccessInfo instance to store \
171  the information about variable accesses.
172  :type var_accesses: \
173  :py:class:`psyclone.core.VariablesAccessInfo`
174 
175  '''
176  cell_ref_name, ref = self.cell_ref_namecell_ref_name(var_accesses)
177  self.psyir_appendpsyir_append(ref)
178  self.appendappend(cell_ref_name)
179 
180  def cell_map(self, var_accesses=None):
181  '''Add cell-map and related cell counts (for inter-grid kernels)
182  to the argument list. If supplied it also stores these accesses to the
183  var_access object.
184 
185  :param var_accesses: optional VariablesAccessInfo instance to store \
186  the information about variable accesses.
187  :type var_accesses: \
188  :py:class:`psyclone.core.VariablesAccessInfo`
189 
190  '''
191  cargs = psyGen.args_filter(self._kern_kern.args, arg_meshes=["gh_coarse"])
192  carg = cargs[0]
193  fargs = psyGen.args_filter(self._kern_kern.args, arg_meshes=["gh_fine"])
194  farg = fargs[0]
195  base_name = "cell_map_" + carg.name
196 
197  # Add the cell map to our argument list
198  cell_ref_name, cell_ref = self.cell_ref_namecell_ref_name(var_accesses)
199  sym = self.append_array_referenceappend_array_reference(base_name, [":", ":", cell_ref],
200  ScalarType.Intrinsic.INTEGER)
201  self.appendappend(f"{sym.name}(:,:,{cell_ref_name})",
202  var_accesses=var_accesses, var_access_name=sym.name)
203 
204  # No. of fine cells per coarse cell in x
205  base_name = f"ncpc_{farg.name}_{carg.name}_x"
206  sym = self.append_integer_referenceappend_integer_reference(base_name)
207  self.appendappend(sym.name, var_accesses)
208  # No. of fine cells per coarse cell in y
209  base_name = f"ncpc_{farg.name}_{carg.name}_y"
210  sym = self.append_integer_referenceappend_integer_reference(base_name)
211  self.appendappend(sym.name, var_accesses)
212  # No. of columns in the fine mesh
213  base_name = f"ncell_{farg.name}"
214  sym = self.append_integer_referenceappend_integer_reference(base_name)
215  self.appendappend(sym.name, var_accesses)
216 
217  def mesh_height(self, var_accesses=None):
218  '''Add mesh height (nlayers) to the argument list and if supplied
219  stores this access in var_accesses.
220 
221  :param var_accesses: optional VariablesAccessInfo instance to store \
222  the information about variable accesses.
223  :type var_accesses: \
224  :py:class:`psyclone.core.VariablesAccessInfo`
225 
226  '''
227  if self._kern_kern.iterates_over not in ["cell_column", "domain"]:
228  return
229  nlayers_symbol = self.append_integer_referenceappend_integer_reference("nlayers")
230  self.appendappend(nlayers_symbol.name, var_accesses)
231  self._nlayers_positions_nlayers_positions.append(self.num_argsnum_args)
232 
233  def scalar(self, scalar_arg, var_accesses=None):
234  '''
235  Add the necessary argument for a scalar quantity as well as an
236  appropriate Symbol to the SymbolTable.
237 
238  :param scalar_arg: the scalar kernel argument.
239  :type scalar_arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
240  :param var_accesses: optional VariablesAccessInfo instance that \
241  stores information about variable accesses.
242  :type var_accesses: \
243  :py:class:`psyclone.core.VariablesAccessInfo`
244 
245  '''
246  super().scalar(scalar_arg, var_accesses)
247  if scalar_arg.is_literal:
248  self.psyir_appendpsyir_append(scalar_arg.psyir_expression())
249  else:
250  sym = self._symtab_symtab.lookup(scalar_arg.name)
251  self.psyir_appendpsyir_append(Reference(sym))
252 
253  # TODO uncomment this method when ensuring we only pass ncell3d once
254  # to any given kernel.
255  # def mesh_ncell3d(self):
256  # ''' Add the number of cells in the full 3D mesh to the argument
257  # list '''
258  # ncell3d_name = self._name_space_manager.create_name(
259  # root_name="ncell_3d", context="PSyVars", label="ncell3d")
260  # self.append(ncell3d_name)
261 
262  def _mesh_ncell2d(self, var_accesses=None):
263  '''Add the number of columns in the mesh to the argument list and if
264  supplied stores this access in var_accesses.
265 
266  :param var_accesses: optional VariablesAccessInfo instance to store \
267  the information about variable accesses.
268  :type var_accesses: \
269  :py:class:`psyclone.core.VariablesAccessInfo`
270 
271  '''
272  sym = self.append_integer_referenceappend_integer_reference("ncell_2d")
273  self.appendappend(sym.name, var_accesses)
274 
275  def _mesh_ncell2d_no_halos(self, var_accesses=None):
276  '''Add the number of columns in the mesh (excluding those in the halo)
277  to the argument list and store this access in var_accesses (if
278  supplied).
279 
280  :param var_accesses: optional VariablesAccessInfo instance to store \
281  the information about variable accesses.
282  :type var_accesses: \
283  :py:class:`psyclone.core.VariablesAccessInfo`
284 
285  '''
286  ncell_symbol = self.append_integer_referenceappend_integer_reference("ncell_2d_no_halos")
287  self.appendappend(ncell_symbol.name, var_accesses)
288 
289  def cma_operator(self, arg, var_accesses=None):
290  '''Add the CMA operator and associated scalars to the argument
291  list and optionally add them to the variable access
292  information.
293 
294  :param arg: the CMA operator argument.
295  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
296  :param var_accesses: optional VariablesAccessInfo instance to store \
297  the information about variable accesses.
298  :type var_accesses: \
299  :py:class:`psyclone.core.VariablesAccessInfo`
300 
301  '''
302  components = ["matrix"]
303  # Avoid circular import:
304  # pylint: disable=import-outside-toplevel
305  from psyclone.dynamo0p3 import DynCMAOperators
306  if arg.function_space_to.orig_name != (arg.function_space_from.
307  orig_name):
308  components += DynCMAOperators.cma_diff_fs_params
309  else:
310  components += DynCMAOperators.cma_same_fs_params
311 
312  const = LFRicConstants()
313  suffix = const.ARG_TYPE_SUFFIX_MAPPING["gh_columnwise_operator"]
314 
315  for component in components:
316  # Matrix takes the access from the declaration of the argument
317  # (i.e. read, write, ...), the rest are always read-only parameters
318  if component == "matrix":
319  # Matrix is a pointer to a 3d array
320  # REAL(KIND=r_solver), pointer:: cma_op1_matrix(:,:,:)
321  # = > null()
322  mode = arg.access
323  sym = self._symtab_symtab.lookup_with_tag(f"{arg.name}:{suffix}")
324  self.psyir_appendpsyir_append(ArrayReference.create(sym, [":", ":", ":"]))
325  else:
326  # All other variables are scalar integers
327  name = self._symtab_symtab.lookup_with_tag(
328  f"{arg.name}:{component}:{suffix}").name
329  mode = AccessType.READ
330  sym = self.append_integer_referenceappend_integer_reference(
331  name, tag=f"{arg.name}:{component}:{suffix}")
332 
333  self.appendappend(sym.name, var_accesses, mode=mode,
334  metadata_posn=arg.metadata_index)
335 
336  def field_vector(self, argvect, var_accesses=None):
337  '''Add the field vector associated with the argument 'argvect' to the
338  argument list. If supplied it also stores these accesses to the
339  var_access object.
340 
341  :param argvect: the field vector to add.
342  :type argvect: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
343  :param var_accesses: optional VariablesAccessInfo instance to store \
344  the information about variable accesses.
345  :type var_accesses: \
346  :py:class:`psyclone.core.VariablesAccessInfo`
347 
348  '''
349  suffix = LFRicConstants().ARG_TYPE_SUFFIX_MAPPING[
350  argvect.argument_type]
351  # The range function below returns values from
352  # 1 to the vector size which is what we
353  # require in our Fortran code
354  for idx in range(1, argvect.vector_size + 1):
355  cmpt_sym = self._symtab_symtab.lookup_with_tag(
356  f"{argvect.name}_{idx}:{suffix}")
357  self.psyir_appendpsyir_append(Reference(cmpt_sym))
358  text = cmpt_sym.name
359  self.appendappend(text, metadata_posn=argvect.metadata_index)
360 
361  if var_accesses is not None:
362  # We add the whole field-vector, not the individual accesses.
363  var_accesses.add_access(Signature(argvect.name), argvect.access,
364  self._kern_kern)
365 
366  def field(self, arg, var_accesses=None):
367  '''Add the field array associated with the argument 'arg' to the
368  argument list. If supplied it also stores this access in var_accesses.
369 
370  :param arg: the field to be added.
371  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
372  :param var_accesses: optional VariablesAccessInfo instance to store
373  the information about variable accesses.
374  :type var_accesses: :py:class:`psyclone.core.VariablesAccessInfo`
375 
376  '''
377  const = LFRicConstants()
378  suffix = const.ARG_TYPE_SUFFIX_MAPPING[arg.argument_type]
379  # Look-up the name of the variable that stores the reference to
380  # the data in this field.
381  sym = self._symtab_symtab.lookup_with_tag(f"{arg.name}:{suffix}")
382  # Add the field data array as being read.
383  self.appendappend(sym.name, var_accesses, var_access_name=sym.name,
384  mode=arg.access, metadata_posn=arg.metadata_index)
385 
386  self.psyir_appendpsyir_append(Reference(sym))
387 
388  def stencil_unknown_extent(self, arg, var_accesses=None):
389  '''Add stencil information to the argument list associated with the
390  argument 'arg' if the extent is unknown. If supplied it also stores
391  this access in var_accesses.
392 
393  :param arg: the kernel argument with which the stencil is associated.
394  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
395  :param var_accesses: optional VariablesAccessInfo instance to store \
396  the information about variable accesses.
397  :type var_accesses: \
398  :py:class:`psyclone.core.VariablesAccessInfo`
399 
400  '''
401  # The extent is not specified in the metadata so pass the value in
402  # Import here to avoid circular dependency
403  # pylint: disable=import-outside-toplevel
404  from psyclone.domain.lfric.lfric_stencils import LFRicStencils
405  var_sym = LFRicStencils.dofmap_size_symbol(self._symtab_symtab, arg)
406  cell_name, cell_ref = self.cell_ref_namecell_ref_name(var_accesses)
407  self.append_array_referenceappend_array_reference(var_sym.name, [cell_ref],
408  ScalarType.Intrinsic.INTEGER,
409  symbol=var_sym)
410  self.appendappend(f"{var_sym.name}({cell_name})", var_accesses,
411  var_access_name=var_sym.name)
412 
413  def stencil_2d_unknown_extent(self, arg, var_accesses=None):
414  '''Add 2D stencil information to the argument list associated with the
415  argument 'arg' if the extent is unknown. If supplied it also stores
416  this access in var_accesses.
417 
418  :param arg: the kernel argument with which the stencil is associated.
419  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
420  :param var_accesses: optional VariablesAccessInfo instance to store \
421  the information about variable accesses.
422  :type var_accesses: \
423  :py:class:`psyclone.core.VariablesAccessInfo`
424 
425  '''
426  # The extent is not specified in the metadata so pass the value in
427  # Import here to avoid circular dependency
428  # pylint: disable=import-outside-toplevel
429  from psyclone.domain.lfric.lfric_stencils import LFRicStencils
430  var_sym = LFRicStencils.dofmap_size_symbol(self._symtab_symtab, arg)
431  cell_name, cell_ref = self.cell_ref_namecell_ref_name(var_accesses)
432  self.append_array_referenceappend_array_reference(var_sym.name, [":", cell_ref],
433  ScalarType.Intrinsic.INTEGER,
434  symbol=var_sym)
435  name = f"{var_sym.name}(:,{cell_name})"
436  self.appendappend(name, var_accesses, var_access_name=var_sym.name)
437 
438  def stencil_2d_max_extent(self, arg, var_accesses=None):
439  '''Add the maximum branch extent for a 2D stencil associated with the
440  argument 'arg' to the argument list. If supplied it also stores this
441  in var_accesses.
442 
443  :param arg: the kernel argument with which the stencil is associated.
444  :type arg: :py:class:`pclone.dynamo0p3.DynKernelArgument`
445  :param var_accesses: optional SingleVariableAccessInfo instance \
446  to store the information about variable accesses.
447  :type var_accesses: \
448  :py:class:`psyclone.core.SingleVariableAccessInfo`
449 
450  '''
451  # The maximum branch extent is not specified in the metadata so pass
452  # the value in.
453  # Import here to avoid circular dependency
454  # pylint: disable=import-outside-toplevel
455  from psyclone.domain.lfric.lfric_stencils import LFRicStencils
456  # TODO #1915, this duplicates code in
457  # LFRicStencils.max_branch_length_name
458  unique_tag = LFRicStencils.stencil_unique_str(arg, "length")
459  root_name = arg.name + "_max_branch_length"
460 
461  sym = self.append_integer_referenceappend_integer_reference(root_name, tag=unique_tag)
462  self.appendappend(sym.name, var_accesses)
463 
464  def stencil_unknown_direction(self, arg, var_accesses=None):
465  '''Add stencil information to the argument list associated with the
466  argument 'arg' if the direction is unknown (i.e. it's being supplied
467  in a variable). If supplied it also stores this access in
468  var_accesses.
469 
470  :param arg: the kernel argument with which the stencil is associated.
471  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
472  :param var_accesses: optional VariablesAccessInfo instance to store \
473  the information about variable accesses.
474  :type var_accesses: \
475  :py:class:`psyclone.core.VariablesAccessInfo`
476 
477  '''
478  # the direction of the stencil is not known so pass the value in
479  name = arg.stencil.direction_arg.varname
480  tag = arg.stencil.direction_arg.text
481  self.append_integer_referenceappend_integer_reference(name, f"AlgArgs_{tag}")
482  self.appendappend(name, var_accesses)
483 
484  def stencil(self, arg, var_accesses=None):
485  '''Add general stencil information associated with the argument 'arg'
486  to the argument list. If supplied it also stores this access in
487  var_accesses.
488 
489  :param arg: the meta-data description of the kernel \
490  argument with which the stencil is associated.
491  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
492  :param var_accesses: optional VariablesAccessInfo instance to store \
493  the information about variable accesses.
494  :type var_accesses: \
495  :py:class:`psyclone.core.VariablesAccessInfo`
496 
497  '''
498  # add in stencil dofmap
499  # Import here to avoid circular dependency
500  # pylint: disable=import-outside-toplevel
501  from psyclone.domain.lfric.lfric_stencils import LFRicStencils
502  var_sym = LFRicStencils.dofmap_symbol(self._symtab_symtab, arg)
503  cell_name, cell_ref = self.cell_ref_namecell_ref_name(var_accesses)
504  self.append_array_referenceappend_array_reference(var_sym.name, [":", ":", cell_ref],
505  ScalarType.Intrinsic.INTEGER,
506  symbol=var_sym)
507  self.appendappend(f"{var_sym.name}(:,:,{cell_name})", var_accesses,
508  var_access_name=var_sym.name)
509 
510  def stencil_2d(self, arg, var_accesses=None):
511  '''Add general 2D stencil information associated with the argument
512  'arg' to the argument list. If supplied it also stores this access in
513  var_accesses.
514 
515  :param arg: the meta-data description of the kernel \
516  argument with which the stencil is associated.
517  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
518  :param var_accesses: optional VariablesAccessInfo instance to store \
519  the information about variable accesses.
520  :type var_accesses: \
521  :py:class:`psyclone.core.VariablesAccessInfo`
522 
523  '''
524  # The stencil_2D differs from the stencil in that the direction
525  # of the branch is baked into the stencil_dofmap array.
526  # The array dimensions are thus (dof_in_cell, cell_in_branch,
527  # branch_in_stencil) where the branch_in_stencil is always ordered
528  # West, South, East, North which is standard in LFRic. This allows
529  # for knowledge of what direction a stencil cell is in relation
530  # to the center even when the stencil is truncated at boundaries.
531  # Import here to avoid circular dependency
532  # pylint: disable=import-outside-toplevel
533  from psyclone.domain.lfric.lfric_stencils import LFRicStencils
534  var_sym = LFRicStencils.dofmap_symbol(self._symtab_symtab, arg)
535  cell_name, cell_ref = self.cell_ref_namecell_ref_name(var_accesses)
536  self.append_array_referenceappend_array_reference(var_sym.name,
537  [":", ":", ":", cell_ref],
538  ScalarType.Intrinsic.INTEGER,
539  symbol=var_sym)
540  name = f"{var_sym.name}(:,:,:,{cell_name})"
541  self.appendappend(name, var_accesses, var_access_name=var_sym.name)
542 
543  def operator(self, arg, var_accesses=None):
544  '''Add the operator arguments to the argument list. If supplied it
545  also stores this access in var_accesses.
546 
547  :param arg: the meta-data description of the operator.
548  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
549  :param var_accesses: optional VariablesAccessInfo instance to store \
550  the information about variable accesses.
551  :type var_accesses: \
552  :py:class:`psyclone.core.VariablesAccessInfo`
553 
554  '''
555  # TODO we should only be including ncell_3d once in the argument
556  # list but this adds it for every operator
557  # This argument is always read only:
558  if arg.data_type == "r_solver_operator_type":
559  op_name = "r_solver_operator"
560  elif arg.data_type == "r_tran_operator_type":
561  op_name = "r_tran_operator"
562  else:
563  op_name = "operator"
564  const = LFRicConstants()
565  operator = const.DATA_TYPE_MAP[op_name]
566  self.append_structure_referenceappend_structure_reference(
567  operator["module"], operator["proxy_type"], ["ncell_3d"],
568  arg.proxy_name_indexed,
569  overwrite_datatype=LFRicTypes("LFRicIntegerScalarDataType")())
570  self.appendappend(arg.proxy_name_indexed + "%ncell_3d", var_accesses,
571  mode=AccessType.READ)
572 
573  sym = self._symtab_symtab.lookup_with_tag(
574  f"{arg.name}:{const.ARG_TYPE_SUFFIX_MAPPING[arg.argument_type]}")
575  self.psyir_appendpsyir_append(Reference(sym))
576  # The access mode of `local_stencil` is taken from the meta-data:
577  self.appendappend(sym.name, var_accesses,
578  mode=arg.access, metadata_posn=arg.metadata_index)
579 
580  def fs_common(self, function_space, var_accesses=None):
581  '''Add function-space related arguments common to LMA operators and
582  fields. If supplied it also stores this access in var_accesses.
583 
584  :param function_space: the function space for which the related \
585  arguments common to LMA operators and fields are added.
586  :type function_space: :py:class:`psyclone.domain.lfric.FunctionSpace`
587  :param var_accesses: optional VariablesAccessInfo instance to store \
588  the information about variable accesses.
589  :type var_accesses: \
590  :py:class:`psyclone.core.VariablesAccessInfo`
591 
592  '''
593  if self._kern_kern.iterates_over not in ["cell_column", "domain"]:
594  return
595  super().fs_common(function_space, var_accesses)
596  self._ndf_positions_ndf_positions.append(
597  KernCallArgList.NdfInfo(position=self.num_argsnum_args,
598  function_space=function_space.orig_name))
599 
600  def fs_compulsory_field(self, function_space, var_accesses=None):
601  '''Add compulsory arguments associated with this function space to
602  the list. If supplied it also stores this access in var_accesses.
603 
604  :param function_space: the function space for which the compulsory \
605  arguments are added.
606  :type function_space: :py:class:`psyclone.domain.lfric.FunctionSpace`
607  :param var_accesses: optional VariablesAccessInfo instance to store \
608  the information about variable accesses.
609  :type var_accesses: \
610  :py:class:`psyclone.core.VariablesAccessInfo`
611 
612  '''
613  sym = self.append_integer_referenceappend_integer_reference(function_space.undf_name)
614  self.appendappend(sym.name, var_accesses)
615 
616  map_name = function_space.map_name
617  if self._kern_kern.iterates_over == 'domain':
618  # This kernel takes responsibility for iterating over cells so
619  # pass the whole dofmap.
620  sym = self.append_array_referenceappend_array_reference(map_name, [":", ":"],
621  ScalarType.Intrinsic.INTEGER)
622  self.appendappend(sym.name, var_accesses, var_access_name=sym.name)
623  else:
624  # Pass the dofmap for the cell column
625  cell_name, cell_ref = self.cell_ref_namecell_ref_name(var_accesses)
626  sym = self.append_array_referenceappend_array_reference(map_name, [":", cell_ref],
627  ScalarType.Intrinsic.INTEGER)
628  self.appendappend(f"{sym.name}(:,{cell_name})",
629  var_accesses, var_access_name=sym.name)
630 
631  def fs_intergrid(self, function_space, var_accesses=None):
632  '''Add function-space related arguments for an intergrid kernel.
633  If supplied it also stores this access in var_accesses.
634 
635  :param function_space: the function space for which to add arguments
636  :type function_space: :py:class:`psyclone.domain.lfric.FunctionSpace`
637  :param var_accesses: optional VariablesAccessInfo instance to store \
638  the information about variable accesses.
639  :type var_accesses: \
640  :py:class:`psyclone.core.VariablesAccessInfo`
641 
642  '''
643  # Is this FS associated with the coarse or fine mesh? (All fields
644  # on a given mesh must be on the same FS.)
645  arg = self._kern_kern.arguments.get_arg_on_space(function_space)
646  if arg.mesh == "gh_fine":
647  # For the fine mesh, we need ndf, undf and the *whole*
648  # dofmap
649  self.fs_commonfs_commonfs_common(function_space, var_accesses=var_accesses)
650  sym = self.append_integer_referenceappend_integer_reference(function_space.undf_name)
651  self.appendappend(sym.name, var_accesses)
652  map_name = function_space.map_name
653  sym = self.append_array_referenceappend_array_reference(map_name, [":", ":"],
654  ScalarType.Intrinsic.INTEGER)
655  self.appendappend(sym.name, var_accesses)
656  else:
657  # For the coarse mesh we only need undf and the dofmap for
658  # the current column
659  self.fs_compulsory_fieldfs_compulsory_fieldfs_compulsory_field(function_space,
660  var_accesses=var_accesses)
661 
662  def basis(self, function_space, var_accesses=None):
663  '''Add basis function information for this function space to the
664  argument list and optionally to the variable access information.
665 
666  :param function_space: the function space for which the basis \
667  function is required.
668  :type function_space: :py:class:`psyclone.domain.lfric.FunctionSpace`
669  :param var_accesses: optional VariablesAccessInfo instance to store \
670  the information about variable accesses.
671  :type var_accesses: \
672  :py:class:`psyclone.core.VariablesAccessInfo`
673 
674  '''
675  for rule in self._kern_kern.qr_rules.values():
676  basis_name = function_space.get_basis_name(qr_var=rule.psy_name)
677  sym = self.append_array_referenceappend_array_reference(basis_name, [":", ":", ":", ":"],
678  ScalarType.Intrinsic.REAL)
679  self.appendappend(sym.name, var_accesses)
680 
681  if "gh_evaluator" in self._kern_kern.eval_shapes:
682  # We are dealing with an evaluator and therefore need as many
683  # basis functions as there are target function spaces.
684  for fs_name in self._kern_kern.eval_targets:
685  # The associated FunctionSpace object is the first item in
686  # the tuple dict entry associated with the name of the target
687  # function space
688  fspace = self._kern_kern.eval_targets[fs_name][0]
689  basis_name = function_space.get_basis_name(on_space=fspace)
690  sym = self.append_array_referenceappend_array_reference(basis_name, [":", ":", ":"],
691  ScalarType.Intrinsic.REAL)
692  self.appendappend(sym.name, var_accesses)
693 
694  def diff_basis(self, function_space, var_accesses=None):
695  '''Add differential basis information for the function space to the
696  argument list. If supplied it also stores this access in
697  var_accesses.
698 
699  :param function_space: the function space for which the differential \
700  basis functions are required.
701  :type function_space: :py:class:`psyclone.domain.lfric.FunctionSpace`
702  :param var_accesses: optional VariablesAccessInfo instance to store \
703  the information about variable accesses.
704  :type var_accesses: \
705  :py:class:`psyclone.core.VariablesAccessInfo`
706 
707  '''
708  for rule in self._kern_kern.qr_rules.values():
709  diff_basis_name = function_space.get_diff_basis_name(
710  qr_var=rule.psy_name)
711  sym = self.append_array_referenceappend_array_reference(diff_basis_name,
712  [":", ":", ":", ":"],
713  ScalarType.Intrinsic.REAL)
714  self.appendappend(sym.name, var_accesses)
715 
716  if "gh_evaluator" in self._kern_kern.eval_shapes:
717  # We are dealing with an evaluator and therefore need as many
718  # basis functions as there are target function spaces.
719  for fs_name in self._kern_kern.eval_targets:
720  # The associated FunctionSpace object is the first item in
721  # the tuple dict entry associated with the name of the target
722  # function space
723  fspace = self._kern_kern.eval_targets[fs_name][0]
724  diff_basis_name = function_space.get_diff_basis_name(
725  on_space=fspace)
726  sym = self.append_array_referenceappend_array_reference(diff_basis_name,
727  [":", ":", ":"],
728  ScalarType.Intrinsic.REAL)
729  self.appendappend(sym.name, var_accesses)
730 
731  def field_bcs_kernel(self, function_space, var_accesses=None):
732  '''Implement the boundary_dofs array fix for a field. If supplied it
733  also stores this access in var_accesses.
734 
735  :param function_space: the function space for which boundary dofs \
736  are required.
737  :type function_space: :py:class:`psyclone.domain.lfric.FunctionSpace`
738  :param var_accesses: optional VariablesAccessInfo instance to store \
739  the information about variable accesses.
740  :type var_accesses: \
741  :py:class:`psyclone.core.VariablesAccessInfo`
742 
743  :raises GenerationError: if the bcs kernel does not contain \
744  a field as argument (but e.g. an operator).
745 
746  '''
747  fspace = None
748  for fspace in self._kern_kern.arguments.unique_fss:
749  if fspace.orig_name == "any_space_1":
750  break
751  farg = self._kern_kern.arguments.get_arg_on_space(fspace)
752  # Sanity check - expect the enforce_bc_code kernel to only have
753  # a field argument.
754  if not farg.is_field:
755  const = LFRicConstants()
756  raise GenerationError(
757  f"Expected an argument of {const.VALID_FIELD_NAMES} type "
758  f"from which to look-up boundary dofs for kernel "
759  f"{self._kern.name} but got '{farg.argument_type}'")
760 
761  base_name = "boundary_dofs_" + farg.name
762  sym = self.append_array_referenceappend_array_reference(base_name, [":", ":"],
763  ScalarType.Intrinsic.INTEGER)
764  self.appendappend(sym.name, var_accesses)
765 
766  def operator_bcs_kernel(self, function_space, var_accesses=None):
767  '''Supply necessary additional arguments for the kernel that
768  applies boundary conditions to a LMA operator. If supplied it
769  also stores this access in var_accesses.
770 
771  :param function_space: unused, only for consistency with base class.
772  :type function_space: :py:class:`psyclone.dynamo3.FunctionSpace`
773  :param var_accesses: optional VariablesAccessInfo instance to store \
774  the information about variable accesses.
775  :type var_accesses: \
776  :py:class:`psyclone.core.VariablesAccessInfo`
777 
778  '''
779  # This kernel has only a single LMA operator as argument.
780  # Checks for this are performed in ArgOrdering.generate()
781  op_arg = self._kern_kern.arguments.args[0]
782  base_name = "boundary_dofs_" + op_arg.name
783  sym = self.append_array_referenceappend_array_reference(base_name, [":", ":"],
784  ScalarType.Intrinsic.INTEGER)
785  self.appendappend(sym.name, var_accesses)
786 
787  def mesh_properties(self, var_accesses=None):
788  '''Provide the kernel arguments required for the mesh properties
789  specified in the kernel metadata. If supplied it also stores this
790  access in var_accesses.
791 
792  :param var_accesses: optional VariablesAccessInfo instance to store \
793  the information about variable accesses.
794  :type var_accesses: \
795  :py:class:`psyclone.core.VariablesAccessInfo`
796 
797  '''
798  if self._kern_kern.mesh.properties:
799  # Avoid circular import:
800  # pylint: disable=import-outside-toplevel
801  from psyclone.dynamo0p3 import LFRicMeshProperties
802  self.extendextend(LFRicMeshProperties(self._kern_kern).
803  kern_args(stub=False, var_accesses=var_accesses,
804  kern_call_arg_list=self))
805 
806  def quad_rule(self, var_accesses=None):
807  '''Add quadrature-related information to the kernel argument list.
808  Adds the necessary arguments to the argument list, and optionally
809  adds variable access information to the var_accesses object.
810 
811  :param var_accesses: optional VariablesAccessInfo instance to store \
812  the information about variable accesses.
813  :type var_accesses: \
814  :py:class:`psyclone.core.VariablesAccessInfo`
815 
816  '''
817  # The QR shapes that this routine supports
818  supported_qr_shapes = ["gh_quadrature_xyoz", "gh_quadrature_edge",
819  "gh_quadrature_face"]
820 
821  for shape, rule in self._kern_kern.qr_rules.items():
822  if shape == "gh_quadrature_xyoz":
823  # XYoZ quadrature requires the number of quadrature points in
824  # the horizontal and in the vertical.
825  self._nqp_positions_nqp_positions.append(
826  {"horizontal": self.num_argsnum_args + 1,
827  "vertical": self.num_argsnum_args + 2})
828  self.extendextend(rule.kernel_args, var_accesses)
829  elif shape == "gh_quadrature_edge":
830  # TODO #705 support transformations supplying the number of
831  # quadrature points for edge quadrature.
832  self.extendextend(rule.kernel_args, var_accesses)
833  elif shape == "gh_quadrature_face":
834  # TODO #705 support transformations supplying the number of
835  # quadrature points for face quadrature.
836  self.extendextend(rule.kernel_args, var_accesses)
837  else:
838  raise NotImplementedError(
839  f"quad_rule: no support implemented for quadrature with a "
840  f"shape of '{shape}'. Supported shapes are: "
841  f"{supported_qr_shapes}.")
842  # Now define the arguments using PSyIR:
843  for arg in rule.kernel_args:
844  # Each rule has a `psy_name` (e.g. qr_xyoz), which is appended
845  # to all variable names (e.g. np_xy_qr_xyoz). Remove this
846  # suffix to get the 'generic' name, from which we derive
847  # the correct type:
848  generic_name = arg[:-len(rule.psy_name)-1]
849  if generic_name in ["np_xy", "np_z", "nfaces", "np_xyz",
850  "nedges"]:
851  # np_xy, np_z, nfaces, np_xyz, nedges are all integers:
852  self.append_integer_referenceappend_integer_reference(arg)
853  elif generic_name in ["weights_xy", "weights_z"]:
854  # 1d arrays:
855  # TODO # 1910: These should be pointers
856  self.append_array_referenceappend_array_reference(arg, [":"],
857  ScalarType.Intrinsic.REAL)
858  elif generic_name in ["weights_xyz"]:
859  # 2d arrays:
860  # TODO #1910: These should be pointers
861  self.append_array_referenceappend_array_reference(arg, [":", ":"],
862  ScalarType.Intrinsic.REAL)
863  else:
864  raise InternalError(f"Found invalid kernel argument "
865  f"'{arg}'.")
866 
867  @property
868  def nlayers_positions(self):
869  ''':returns: the position(s) in the argument list of the \
870  variable(s) that passes the number of layers. The generate \
871  method must be called first.
872  :rtype: list of int.
873 
874  :raises InternalError: if the generate() method has not been called.
875 
876  '''
877  if not self._generate_called_generate_called:
878  raise InternalError(
879  "KernCallArgList: the generate() method should be called "
880  "before the nlayers_positions() method")
881  return self._nlayers_positions_nlayers_positions
882 
883  @property
884  def nqp_positions(self):
885  ''':return: the positions in the argument list of the variables that \
886  pass the number of quadrature points. The number and type of \
887  these will change depending on the type of quadrature. A list \
888  of dictionaries is returned with the quadrature types \
889  being the keys to the dictionaries and their position in the \
890  argument list being the values. At the moment only XYoZ is \
891  supported (which has horizontal and vertical quadrature \
892  points). The generate method must be called first.
893  :rtype: [{str: int, ...}]
894 
895  :raises InternalError: if the generate() method has not been \
896  called.
897 
898  '''
899  if not self._generate_called_generate_called:
900  raise InternalError(
901  "KernCallArgList: the generate() method should be called "
902  "before the nqp_positions() method")
903  return self._nqp_positions_nqp_positions
904 
905  @property
906  def ndf_positions(self):
907  ''':return: the position(s) in the argument list and the function \
908  space(s) associated with the variable(s) that pass(es) the \
909  number of degrees of freedom for the function space. The \
910  generate method must be called first.
911  :rtype: list of namedtuple (position=int, function_space=str).
912 
913  :raises InternalError: if the generate() method has not been \
914  called.
915 
916  '''
917  if not self._generate_called_generate_called:
918  raise InternalError(
919  "KernCallArgList: the generate() method should be called "
920  "before the ndf_positions() method")
921  return self._ndf_positions_ndf_positions
922 
923  def cell_ref_name(self, var_accesses=None):
924  '''Utility routine which determines whether to return the cell value
925  or the colourmap lookup value. If supplied it also stores this access
926  in var_accesses.
927 
928  :param var_accesses: optional VariablesAccessInfo instance to store \
929  the information about variable accesses.
930  :type var_accesses: \
931  :py:class:`psyclone.core.VariablesAccessInfo`
932 
933  :returns: the Fortran code needed to access the current cell index.
934  :rtype: Tuple[str, py:class:`psyclone.psyir.nodes.Reference`]
935 
936  '''
937  cell_sym = self._symtab_symtab.find_or_create_integer_symbol(
938  "cell", tag="cell_loop_idx")
939  if self._kern_kern.is_coloured():
940  colour_sym = self._symtab_symtab.find_or_create_integer_symbol(
941  "colour", tag="colours_loop_idx")
942  if self._kern_kern.is_intergrid:
943  tag = None
944  else:
945  # If there is only one colourmap we need to specify the tag
946  # to make sure we get the right symbol.
947  tag = "cmap"
948  array_ref = self.get_array_referenceget_array_reference(self._kern_kern.colourmap,
949  [Reference(colour_sym),
950  Reference(cell_sym)],
951  ScalarType.Intrinsic.INTEGER,
952  tag=tag)
953  if var_accesses is not None:
954  var_accesses.add_access(Signature(colour_sym.name),
955  AccessType.READ, self._kern_kern)
956  var_accesses.add_access(Signature(cell_sym.name),
957  AccessType.READ, self._kern_kern)
958  var_accesses.add_access(Signature(array_ref.name),
959  AccessType.READ,
960  self._kern_kern, ["colour", "cell"])
961 
962  return (self._kern_kern.colourmap + "(colour,cell)",
963  array_ref)
964 
965  if var_accesses is not None:
966  var_accesses.add_access(Signature("cell"), AccessType.READ,
967  self._kern_kern)
968 
969  return (cell_sym.name, Reference(cell_sym))
970 
971 
972 # ============================================================================
973 # For automatic documentation creation:
974 __all__ = ["KernCallArgList"]
def append_array_reference(self, array_name, indices, intrinsic_type, tag=None, symbol=None)
def append(self, var_name, var_accesses=None, var_access_name=None, mode=AccessType.READ, metadata_posn=None)
def fs_common(self, function_space, var_accesses=None)
def fs_compulsory_field(self, function_space, var_accesses=None)
def append_integer_reference(self, name, tag=None)
def extend(self, list_var_name, var_accesses=None, mode=AccessType.READ, list_metadata_posn=None)
def get_array_reference(self, array_name, indices, intrinsic_type, tag=None, symbol=None)
def field_bcs_kernel(self, function_space, var_accesses=None)
def get_user_type(self, module_name, user_type, name, tag=None)
def scalar(self, scalar_arg, var_accesses=None)
def fs_common(self, function_space, var_accesses=None)
def operator_bcs_kernel(self, function_space, var_accesses=None)
def append_structure_reference(self, module_name, user_type, member_list, name, tag=None, overwrite_datatype=None)
def basis(self, function_space, var_accesses=None)
def diff_basis(self, function_space, var_accesses=None)
def fs_intergrid(self, function_space, var_accesses=None)
def fs_compulsory_field(self, function_space, var_accesses=None)