Reference Guide  2.5.0
kern_call_acc_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 all of the data references
39 that must be copied over to a GPU before executing the kernel. Ordering
40 of the parameters does not matter apart from where we have members of
41 derived types. In that case, the derived type itself must be specified
42 first before any members.
43 '''
44 
45 from psyclone import psyGen
46 from psyclone.domain.lfric import KernCallArgList
47 from psyclone.errors import InternalError
48 
49 
51  # TODO: #845 Check that all implicit variables have the right type.
52  '''Kernel call arguments that need to be declared by OpenACC
53  directives. KernCallArgList only needs to be specialised
54  where modified, or additional, arguments are required.
55  Scalars are apparently not required but it is valid in
56  OpenACC to include them and requires less specialisation
57  to keep them in.
58 
59  '''
60  def cell_map(self, var_accesses=None):
61  '''Add cell-map to the list of required arrays.
62 
63  :param var_accesses: optional VariablesAccessInfo instance to store
64  the information about variable accesses.
65  :type var_accesses: Optional[
66  :py:class:`psyclone.core.VariablesAccessInfo`]
67 
68  '''
69  cargs = psyGen.args_filter(self._kern_kern.args, arg_meshes=["gh_coarse"])
70  if len(cargs) > 1:
71  raise InternalError(
72  f"An LFRic intergrid kernel should have only one coarse mesh "
73  f"but '{self._kern.name}' has {len(cargs)}")
74  carg = cargs[0]
75  # Add the cell map to our argument list
76  base_name = "cell_map_" + carg.name
77  self.appendappend(base_name)
78  # We'll need the current cell to index into this cell map.
79  self.cell_positioncell_positioncell_positioncell_position(var_accesses)
80 
81  def cell_position(self, var_accesses=None):
82  '''Adds a cell argument to the argument list and if supplied stores
83  this access in var_accesses. Although normally just a scalar, the cell
84  argument may actually require a lookup from a colour map array. Either
85  way, this method adds the name of the variable to the argument list.
86 
87  :param var_accesses: optional VariablesAccessInfo instance to store
88  the information about variable accesses.
89  :type var_accesses: Optional[
90  :py:class:`psyclone.core.VariablesAccessInfo`]
91 
92  '''
93  _, ref = self.cell_ref_namecell_ref_name(var_accesses)
94  self.appendappend(ref.symbol.name)
95 
96  def stencil(self, arg, var_accesses=None):
97  '''Add general stencil information associated with the argument 'arg'
98  to the argument list. OpenACC requires the full dofmap to be
99  specified. If supplied it also stores this access in var_accesses.
100 
101  :param arg: the meta-data description of the kernel argument with
102  which the stencil is associated.
103  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
104  :param var_accesses: optional VariablesAccessInfo instance to store
105  the information about variable accesses.
106  :type var_accesses: Optional[
107  :py:class:`psyclone.core.VariablesAccessInfo`]
108 
109  '''
110  # Import here to avoid circular dependency
111  # pylint: disable=import-outside-toplevel
112  from psyclone.domain.lfric.lfric_stencils import LFRicStencils
113  var_name = LFRicStencils.dofmap_symbol(self._kern_kern.root.symbol_table,
114  arg).name
115  self.appendappend(var_name, var_accesses)
116 
117  def stencil_2d(self, arg, var_accesses=None):
118  '''Add general 2D stencil information associated with the argument
119  'arg' to the argument list. OpenACC requires the full dofmap to be
120  specified. If supplied it also stores this access in var_accesses.This
121  method passes through to the stencil method.
122 
123  :param arg: the meta-data description of the kernel
124  argument with which the stencil is associated.
125  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
126  :param var_accesses: optional VariablesAccessInfo instance to store
127  the information about variable accesses.
128  :type var_accesses: Optional[
129  :py:class:`psyclone.core.VariablesAccessInfo`]
130 
131  '''
132  self.stencilstencilstencilstencil(arg, var_accesses)
133 
134  def stencil_unknown_extent(self, arg, var_accesses=None):
135  '''Add stencil information to the argument list associated with the
136  argument 'arg' if the extent is unknown. If supplied it also stores
137  this access in var_accesses.
138 
139  :param arg: the kernel argument with which the stencil is associated.
140  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
141  :param var_accesses: optional VariablesAccessInfo instance to store
142  the information about variable accesses.
143  :type var_accesses: Optional[
144  :py:class:`psyclone.core.VariablesAccessInfo`]
145 
146  '''
147  # The extent is not specified in the metadata so pass the value in
148  # Import here to avoid circular dependency
149  # pylint: disable=import-outside-toplevel
150  from psyclone.domain.lfric.lfric_stencils import LFRicStencils
151  name = LFRicStencils.dofmap_size_symbol(self._kern_kern.root.symbol_table,
152  arg).name
153  self.appendappend(name, var_accesses)
154 
155  def stencil_2d_unknown_extent(self, arg, var_accesses=None):
156  '''Add 2D stencil information to the argument list associated with the
157  argument 'arg' if the extent is unknown. If supplied it also stores
158  this access in var_accesses. This method passes through to the
159  stencil_unknown_extent method.
160 
161  :param arg: the kernel argument with which the stencil is associated.
162  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
163  :param var_accesses: optional VariablesAccessInfo instance to store
164  the information about variable accesses.
165  :type var_accesses: Optional[
166  :py:class:`psyclone.core.VariablesAccessInfo`]
167 
168  '''
169  self.stencil_unknown_extentstencil_unknown_extentstencil_unknown_extentstencil_unknown_extent(arg, var_accesses)
170 
171  def operator(self, arg, var_accesses=None):
172  '''Add the operator arguments if they have not already been
173  added. OpenACC requires the derived type and the dereferenced
174  data to be specified. If supplied it also stores this access in
175  var_accesses.
176 
177  :param arg: the meta-data description of the operator.
178  :type arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
179  :param var_accesses: optional VariablesAccessInfo instance to store
180  the information about variable accesses.
181  :type var_accesses: Optional[
182  :py:class:`psyclone.core.VariablesAccessInfo`]
183 
184  '''
185  # In case of OpenACC we do not want to transfer the same
186  # data to GPU twice.
187  if arg.proxy_name_indexed not in self.arglistarglist:
188  self.appendappend(arg.proxy_name_indexed, var_accesses)
189  # This adds ncell_3d and local_stencil after the derived type:
190  super().operator(arg, var_accesses)
191 
192  def fs_compulsory_field(self, function_space, var_accesses=None):
193  '''Add compulsory arguments associated with this function space to
194  the list. OpenACC requires the full function-space map
195  to be specified. If supplied it also stores this access in
196  var_accesses.
197 
198  :param function_space: the function space for which the compulsory
199  arguments are added.
200  :type function_space: :py:class:`psyclone.domain.lfric.FunctionSpace`
201  :param var_accesses: optional VariablesAccessInfo instance to store
202  the information about variable accesses.
203  :type var_accesses: Optional[
204  :py:class:`psyclone.core.VariablesAccessInfo`]
205 
206  '''
207  if self._kern_kern.iterates_over != "cell_column":
208  return
209  self.appendappend(function_space.undf_name, var_accesses)
210  # The base class only adds one dimension to the list, while OpenACC
211  # needs the whole field, so we cannot call the base class.
212  self.appendappend(function_space.map_name, var_accesses)
213 
214  def fs_intergrid(self, function_space, var_accesses=None):
215  '''Add arrays that need to be uploaded for inter-grid kernels.
216  These arrays contain the mapping between fine and coarse meshes.
217 
218  :param function_space: the function space associated with the mesh.
219  :type function_space: :py:class:`psyclone.domain.lfric.FunctionSpace`
220  :param var_accesses: optional VariablesAccessInfo instance to store
221  the information about variable accesses.
222  :type var_accesses: Optional[
223  :py:class:`psyclone.core.VariablesAccessInfo`]
224 
225  '''
226  # Is this FS associated with the coarse or fine mesh? (All fields
227  # on a given mesh must be on the same FS.)
228  arg = self._kern_kern.arguments.get_arg_on_space(function_space)
229  if arg.mesh == "gh_fine":
230  # For the fine mesh, we need the *whole* dofmap
231  map_name = function_space.map_name
232  self.appendappend(map_name, var_accesses)
233  else:
234  # For the coarse mesh we only need undf and the dofmap for
235  # the current column
236  self.fs_compulsory_fieldfs_compulsory_fieldfs_compulsory_fieldfs_compulsory_field(function_space,
237  var_accesses=var_accesses)
238 
239  def scalar(self, scalar_arg, var_accesses=None):
240  '''
241  Override the default implementation as there's no need to specify
242  scalars for an OpenACC data region.
243 
244  :param scalar_arg: the kernel argument.
245  :type scalar_arg: :py:class:`psyclone.dynamo0p3.DynKernelArgument`
246  :param var_accesses: optional VariablesAccessInfo instance that
247  stores information about variable accesses.
248  :type var_accesses: Optional[
249  :py:class:`psyclone.core.VariablesAccessInfo`]
250 
251  '''
252 
253 
254 # ============================================================================
255 # For automatic documentation creation:
256 __all__ = ["KernCallAccArgList"]
def append(self, var_name, var_accesses=None, var_access_name=None, mode=AccessType.READ, metadata_posn=None)
def cell_position(self, var_accesses=None)
def fs_compulsory_field(self, function_space, var_accesses=None)
def stencil(self, arg, var_accesses=None)
def stencil_unknown_extent(self, arg, var_accesses=None)
def fs_compulsory_field(self, function_space, var_accesses=None)
def fs_compulsory_field(self, function_space, var_accesses=None)