File: MPI_Scan.3.rst

package info (click to toggle)
openmpi 5.0.8-10
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 201,692 kB
  • sloc: ansic: 613,078; makefile: 42,351; sh: 11,194; javascript: 9,244; f90: 7,052; java: 6,404; perl: 5,179; python: 1,859; lex: 740; fortran: 61; cpp: 20; tcl: 12
file content (246 lines) | stat: -rw-r--r-- 6,981 bytes parent folder | download | duplicates (4)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
.. _mpi_scan:


MPI_Scan
========

.. include_body

:ref:`MPI_Scan`, :ref:`MPI_Iscan`, :ref:`MPI_Scan_init` - Computes an inclusive scan
(partial reduction)


SYNTAX
------


C Syntax
^^^^^^^^

.. code-block:: c

   #include <mpi.h>

   int MPI_Scan(const void *sendbuf, void *recvbuf, int count,
                MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)

   int MPI_Iscan(const void *sendbuf, void *recvbuf, int count,
                 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm,
                 MPI_Request *request)

   int MPI_Scan_init(const void *sendbuf, void *recvbuf, int count,
                 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm,
                 MPI_Info info, MPI_Request *request)


Fortran Syntax
^^^^^^^^^^^^^^

.. code-block:: fortran

   USE MPI
   ! or the older form: INCLUDE 'mpif.h'
   MPI_SCAN(SENDBUF, RECVBUF, COUNT, DATATYPE, OP, COMM, IERROR)
   	<type>	SENDBUF(*), RECVBUF(*)
   	INTEGER	COUNT, DATATYPE, OP, COMM, IERROR

   MPI_ISCAN(SENDBUF, RECVBUF, COUNT, DATATYPE, OP, COMM, REQUEST, IERROR)
   	<type>	SENDBUF(*), RECVBUF(*)
   	INTEGER	COUNT, DATATYPE, OP, COMM, REQUEST, IERROR

   MPI_SCAN_INIT(SENDBUF, RECVBUF, COUNT, DATATYPE, OP, COMM, INFO, REQUEST, IERROR)
   	<type>	SENDBUF(*), RECVBUF(*)
   	INTEGER	COUNT, DATATYPE, OP, COMM, INFO, REQUEST, IERROR


Fortran 2008 Syntax
^^^^^^^^^^^^^^^^^^^

.. code-block:: fortran

   USE mpi_f08
   MPI_Scan(sendbuf, recvbuf, count, datatype, op, comm, ierror)
   	TYPE(*), DIMENSION(..), INTENT(IN) :: sendbuf
   	TYPE(*), DIMENSION(..) :: recvbuf
   	INTEGER, INTENT(IN) :: count
   	TYPE(MPI_Datatype), INTENT(IN) :: datatype
   	TYPE(MPI_Op), INTENT(IN) :: op
   	TYPE(MPI_Comm), INTENT(IN) :: comm
   	INTEGER, OPTIONAL, INTENT(OUT) :: ierror

   MPI_Iscan(sendbuf, recvbuf, count, datatype, op, comm, request, ierror)
   	TYPE(*), DIMENSION(..), INTENT(IN), ASYNCHRONOUS :: sendbuf
   	TYPE(*), DIMENSION(..), ASYNCHRONOUS :: recvbuf
   	INTEGER, INTENT(IN) :: count
   	TYPE(MPI_Datatype), INTENT(IN) :: datatype
   	TYPE(MPI_Op), INTENT(IN) :: op
   	TYPE(MPI_Comm), INTENT(IN) :: comm
   	TYPE(MPI_Request), INTENT(OUT) :: request
   	INTEGER, OPTIONAL, INTENT(OUT) :: ierror

   MPI_Scan_init(sendbuf, recvbuf, count, datatype, op, comm, info, request, ierror)
   	TYPE(*), DIMENSION(..), INTENT(IN), ASYNCHRONOUS :: sendbuf
   	TYPE(*), DIMENSION(..), ASYNCHRONOUS :: recvbuf
   	INTEGER, INTENT(IN) :: count
   	TYPE(MPI_Datatype), INTENT(IN) :: datatype
   	TYPE(MPI_Op), INTENT(IN) :: op
   	TYPE(MPI_Comm), INTENT(IN) :: comm
   	TYPE(MPI_Info), INTENT(IN) :: info
   	TYPE(MPI_Request), INTENT(OUT) :: request
   	INTEGER, OPTIONAL, INTENT(OUT) :: ierror


INPUT PARAMETERS
----------------
* ``sendbuf``: Send buffer (choice).
* ``count``: Number of elements in input buffer (integer).
* ``datatype``: Data type of elements of input buffer (handle).
* ``op``: Operation (handle).
* ``comm``: Communicator (handle).
* ``info``: Info (handle, persistent only)

OUTPUT PARAMETERS
-----------------
* ``recvbuf``: Receive buffer (choice).
* ``request``: Request (handle, non-blocking only).
* ``ierror``: Fortran only: Error status (integer).

DESCRIPTION
-----------

:ref:`MPI_Scan` is used to perform an inclusive prefix reduction on data
distributed across the calling processes. The operation returns, in the
*recvbuf* of the process with rank i, the reduction (calculated
according to the function *op*) of the values in the *sendbuf*\ s of
processes with ranks 0, ..., i (inclusive). The type of operations
supported, their semantics, and the constraints on send and receive
buffers are as for :ref:`MPI_Reduce`.


EXAMPLE
-------

This example uses a user-defined operation to produce a segmented scan.
A segmented scan takes, as input, a set of values and a set of logicals,
where the logicals delineate the various segments of the scan. For
example,

::

   values     v1      v2      v3      v4      v5      v6      v7      v8
   logicals   0       0       1       1       1       0       0       1
   result     v1    v1+v2     v3    v3+v4  v3+v4+v5   v6    v6+v7     v8

The result for rank j is thus the sum v(i) + ... + v(j), where i is the
lowest rank such that for all ranks n, i <= n <= j, logical(n) =
logical(j). The operator that produces this effect is

::

         [ u ]     [ v ]     [ w ]
         [   ]  o  [   ]  =  [   ]
         [ i ]     [ j ]     [ j ]

   where

       ( u + v if i = j w = ( ( v if i != j

Note that this is a noncommutative operator. C code that implements it
is given below.

.. code-block:: c

   	typedef struct {
   		double val;
   		int log;
   	} SegScanPair;

   	/*
   	 * the user-defined function
   	 */
   	void segScan(SegScanPair *in, SegScanPair *inout, int *len,
   		MPI_Datatype *dptr)
   	{
   		int i;
   		SegScanPair c;

   		for (i = 0; i < *len; ++i) {
   			if (in->log == inout->log)
   				c.val = in->val + inout->val;
   			else
   				c.val = inout->val;

   			c.log = inout->log;
   			*inout = c;
   			in++;
   			inout++;
   		}
   	}

Note that the inout argument to the user-defined function corresponds to
the right-hand operand of the operator. When using this operator, we
must be careful to specify that it is noncommutative, as in the
following:

.. code-block:: c

   	int			i, base;
   	SeqScanPair	a, answer;
   	MPI_Op		myOp;
   	MPI_Datatype	type[2] = {MPI_DOUBLE, MPI_INT};
   	MPI_Aint		disp[2];
   	int			blocklen[2] = {1, 1};
   	MPI_Datatype	sspair;

   	/*
   	 * explain to MPI how type SegScanPair is defined
   	 */
   	MPI_Get_address(a, disp);
   	MPI_Get_address(a.log, disp + 1);
   	base = disp[0];
   	for (i = 0; i < 2; ++i)
   		disp[i] -= base;
   	MPI_Type_struct(2, blocklen, disp, type, &sspair);
   	MPI_Type_commit(&sspair);

   	/*
   	 * create the segmented-scan user-op
   	 * noncommutative - set commute (arg 2) to 0
   	 */
   	MPI_Op_create((MPI_User_function *)segScan, 0, &myOp);
   	...
   	MPI_Scan(a, answer, 1, sspair, myOp, comm);


USE OF IN-PLACE OPTION
----------------------

When the communicator is an intracommunicator, you can perform a
scanning operation in place (the output buffer is used as the input
buffer). Use the variable MPI_IN_PLACE as the value of the *sendbuf*
argument. The input data is taken from the receive buffer and replaced
by the output data.


NOTES ON COLLECTIVE OPERATIONS
------------------------------

The reduction functions of type MPI_Op do not return an error value. As
a result, if the functions detect an error, all they can do is either
call :ref:`MPI_Abort` or silently skip the problem. Thus, if the error handler
is changed from MPI_ERRORS_ARE_FATAL to something else (e.g.,
MPI_ERRORS_RETURN), then no error may be indicated.

The reason for this is the performance problems in ensuring that all
collective routines return the same error value.


ERRORS
------

.. include:: ./ERRORS.rst

.. seealso::
   * :ref:`MPI_Exscan`
   * :ref:`MPI_Op_create`
   * :ref:`MPI_Reduce`