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
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
|
/*@@
@file ReductionMax.c
@date
@author Tom Goodale, Yaakoub Y El Khamra
@desc
Defines the reduction operator to get the average
of an arbitrary array.
@enddesc
@version $Id$
@@*/
#include <stdlib.h>
#include <string.h>
#include "local_reductions.h"
#include "Max_Functions.h"
static const char *rcsid = "$Id$";
CCTK_FILEVERSION(CactusBase_LocalReduce_ReductionMax_c);
/* Define the reduction operations */
/* local function prototypes */
static int ReductionMax (int N_dims, int operator_handle,
int param_table_handle, int N_input_arrays,
const CCTK_INT input_array_dims[],
const CCTK_INT input_array_type_codes[],
const void *const input_arrays[],
int M_output_numbers,
const CCTK_INT output_number_type_codes[],
void * const output_numbers[]);
/*@@
@routine LocalReduce_Max
@author Tom Goodale, Yaakoub Y El Khamra
@date
@desc
@enddesc
@history
@endhistory
@var N_dims
@vdesc number of dimensions in the *reduction*
@vtype int
@vio in
@endvar
@var operator_handle
@vdesc operator handle specificies the type of reduction we will perform
@vtype int
@vio in
@endvar
@var param_table_handle
@vdesc handle to "parameter table", a key-value table
@vtype int
@vio in
@endvar
@var N_input_arrays
@vdesc number of input arrays
@vtype int
@vio in
@endvar
@var input_array_dims
@vdesc array of input array dimensions (common to all input arrays)
@vtype const CCTK_INT
@vio in
@endvar
@var input_array_type_codes
@vdesc array of CCTK_VARIABLE_* codes giving data types of input arrays
@vtype const CCTK_INT
@vio in
@endvar
@var input_arrays
@vdesc array of pointers to input arrays
@vtype const void *const
@vio in
@endvar
@var M_output_numbers
@vdesc
@vtype int
@vio in
@endvar
@var output_number_type_codes
@vdesc array of CCTK_VARIABLE_* codes giving data types of output numbers
@vtype const CCTK_INT
@vio in
@endvar
@var output_numbers
@vdesc array[M_output_numbers] of pointers to output numbers[M_reduce_numbers]
@vtype void *const
@vio in
@endvar
@@*/
int LocalReduce_Max (int N_dims, int operator_handle,
int param_table_handle, int N_input_arrays,
const CCTK_INT input_array_dims[],
const CCTK_INT input_array_type_codes[],
const void *const input_arrays[],
int M_output_numbers,
const CCTK_INT output_number_type_codes[],
void * const output_numbers[])
{
return (LocalReduce_Reduce (N_dims, operator_handle,
param_table_handle, N_input_arrays,
input_array_dims, input_array_type_codes,
input_arrays, M_output_numbers,
output_number_type_codes, output_numbers,
ReductionMax));
}
/*****************************************************************************/
/* local functions */
/*****************************************************************************/
/*@@
@routine ReductionMax
@date
@author Tom Goodale, Yaakoub Y El Khamra
@desc Returns the average of a distributed array with
'num_points' elements. Global reduction is done element-wise
(num_outvals == 1) or on the results of the local reductions.
@enddesc
@@*/
static int ReductionMax (int N_dims, int operator_handle,
int param_table_handle, int N_input_arrays,
const CCTK_INT input_array_dims[],
const CCTK_INT input_array_type_codes[],
const void *const input_arrays[],
int M_output_numbers,
const CCTK_INT output_number_type_codes[],
void * const output_numbers[])
{
/* utility variables */
int i, j, num_points;
int ierr;
int * iters_per_dim;
int global_calling = 0;
/* indices to hold the temp indices of size N_dims and iteration indices*/
int * indices;
int * actual_indices;
int * actual_iters_per_dim;
int max_iter = 0;
int max_index = 1;
/* data pointer offset and strides declared here */
CCTK_INT * input_array_offsets;
CCTK_INT * input_array_strides;
CCTK_INT * input_array_min_subscripts;
CCTK_INT * input_array_max_subscripts;
/* weight variables declared here */
int weight_on = 0; /* weight is by default off=0 */
void const * weight; /* pointer to the weight variable */
CCTK_REAL weight_sum;
/* prevent warnings for unused vars */
(void)(operator_handle + 0);
/* set the number of points */
num_points = 0;
/* allocate memory for iters_per_dim */
iters_per_dim = (int *)malloc(N_dims * sizeof(int));
/* allocate then initialize the values of the strides and subscripts */
indices = (int *)malloc (N_dims * sizeof(int));
actual_indices = (int *)malloc (N_dims * sizeof(int));
actual_iters_per_dim = (int *)malloc (N_dims * sizeof(int));
/* allocate then initialize the values of the strides and subscripts */
input_array_offsets = (CCTK_INT *)malloc (N_input_arrays * sizeof(CCTK_INT));
input_array_strides = (CCTK_INT *) malloc (N_dims * sizeof(CCTK_INT));
input_array_min_subscripts = (CCTK_INT *) malloc (N_dims * sizeof(CCTK_INT));
input_array_max_subscripts = (CCTK_INT *) malloc (N_dims * sizeof(CCTK_INT));
for (i = 0; i<N_input_arrays; i++)
{
input_array_offsets[i] = 0;
}
for (i = 0; i<N_dims; i++)
{
input_array_strides[i] = 1;
input_array_min_subscripts[i] = 0;
input_array_max_subscripts[i] = input_array_dims[i];
max_index *= input_array_max_subscripts[i];
}
/* for strides and subscripts get values from param table (it they exist) */
if ( Util_TableQueryNKeys(param_table_handle) != 0)
{
ierr = Util_TableGetGenericArray(param_table_handle, CCTK_VARIABLE_INT,
N_input_arrays, input_array_offsets, "input_array_offsets");
ierr = Util_TableGetGenericArray(param_table_handle, CCTK_VARIABLE_INT,
N_dims, input_array_strides, "input_array_strides");
ierr = Util_TableGetGenericArray(param_table_handle, CCTK_VARIABLE_INT,
N_dims, input_array_min_subscripts, "input_array_min_subscripts");
ierr = Util_TableGetGenericArray(param_table_handle, CCTK_VARIABLE_INT,
N_dims, input_array_max_subscripts, "input_array_max_subscripts");
}
/* for weights get values from param table (it they exist) */
if ( Util_TableQueryNKeys(param_table_handle) != 0)
{
ierr = 0;
/* Get the weight_on parameter and the weight_var_index parameter */
ierr = Util_TableGetInt(param_table_handle, &weight_on, "weight_on");
ierr = Util_TableGetPointerToConst(param_table_handle, &weight, "weight");
/* Need to add more checking here for size, storage */
}
/* reduction maps an array to a single value of the same type */
if (M_output_numbers != N_input_arrays)
{
CCTK_WARN (1, "Average reduction returns a single value\n \
for each input array\n");
return (-1);
}
/* set the indices to their minimum values */
max_iter = 1;
for (j = 0; j <N_dims; j++)
{
indices [j] = 0;
actual_indices[j] = input_array_min_subscripts[j];
actual_iters_per_dim [j] = (int) (input_array_max_subscripts[j] - input_array_min_subscripts[j]);
iters_per_dim [j] = (int) ((input_array_max_subscripts[j] - input_array_min_subscripts[j])/input_array_strides[j]);
max_iter *= iters_per_dim [j];
}
for (i = 0; i < N_input_arrays; i++)
{
/* Do the type matching */
switch (input_array_type_codes[i])
{
/* in values type switches*/
case CCTK_VARIABLE_BYTE:
ierr = LocalReduce_Max_BYTE(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
case CCTK_VARIABLE_INT:
ierr = LocalReduce_Max_INT(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#ifdef CCTK_INT1
case CCTK_VARIABLE_INT1:
ierr = LocalReduce_Max_INT1(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#endif
#ifdef CCTK_INT2
case CCTK_VARIABLE_INT2:
ierr = LocalReduce_Max_INT2(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#endif
#ifdef CCTK_INT4
case CCTK_VARIABLE_INT4:
ierr = LocalReduce_Max_INT4(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#endif
#ifdef CCTK_INT8
case CCTK_VARIABLE_INT8:
ierr = LocalReduce_Max_INT8(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#endif
case CCTK_VARIABLE_REAL:
ierr = LocalReduce_Max_REAL(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#ifdef CCTK_REAL4
case CCTK_VARIABLE_REAL4:
ierr = LocalReduce_Max_REAL4(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#endif
#ifdef CCTK_REAL8
case CCTK_VARIABLE_REAL8:
ierr = LocalReduce_Max_REAL8(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#endif
#ifdef CCTK_REAL16
case CCTK_VARIABLE_REAL16:
ierr = LocalReduce_Max_REAL16(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#endif
case CCTK_VARIABLE_COMPLEX:
ierr = LocalReduce_Max_COMPLEX(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#ifdef CCTK_COMPLEX8
case CCTK_VARIABLE_COMPLEX8:
ierr = LocalReduce_Max_COMPLEX8(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#endif
#ifdef CCTK_COMPLEX16
case CCTK_VARIABLE_COMPLEX16:
ierr = LocalReduce_Max_COMPLEX16(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#endif
#ifdef CCTK_COMPLEX32
case CCTK_VARIABLE_COMPLEX32:
ierr = LocalReduce_Max_COMPLEX32(i, weight_on, weight, input_array_offsets, indices, max_iter, actual_indices, input_array_strides, input_array_min_subscripts, input_array_dims, num_points, actual_iters_per_dim, iters_per_dim, N_dims, input_arrays, output_number_type_codes, output_numbers, param_table_handle);
break;
#endif
}
}
/* Get the values of num_points and weight_sum */
ierr = Util_TableGetInt(param_table_handle, &num_points, "num_points");
num_points--;
/* store the number of points in the paramater table and perform division */
ierr = Util_TableGetInt(param_table_handle, &global_calling, "global_calling");
if ( global_calling != 0)
{
ierr = Util_TableSetInt(param_table_handle, num_points, "num_points");
ierr = Util_TableSetInt(param_table_handle, 1,"global_operation");
}
/* free memory */
free (iters_per_dim);
free (indices);
free (actual_indices);
free (actual_iters_per_dim);
free (input_array_offsets);
free (input_array_strides);
free (input_array_min_subscripts);
free (input_array_max_subscripts);
return (0);
}
|