summaryrefslogtreecommitdiff
path: root/src/comm/Reduction.c
blob: bef214b520d5c946a1462de5849e9d047025da73 (plain)
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
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
 /*@@
   @file      Reduction.c
   @date      1999/05/13
   @author    Gabrielle Allen
   @desc
              This file contains routines to deal with registering and
              using functions providing reduction operations.
   @enddesc
   @version   $Id$
 @@*/

/* #define DEBUG_REDUCTION 1 */

#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <stdarg.h>

#include "cctk_Flesh.h"
#include "cctk_FortranString.h"
#include "cctk_Groups.h"
#include "cctk_Reduction.h"
#include "cctk_WarnLevel.h"
#include "cctk_ActiveThorns.h"

#include "StoreHandledData.h"

static const char *rcsid = "$Header$";

CCTK_FILEVERSION(comm_Reduction_c);

/********************************************************************
 ********************    External Routines   ************************
 ********************************************************************/
/* prototypes for external C routines are declared in header cctk_Reduce.h
   here only follow the fortran wrapper prototypes */

void CCTK_FCALL CCTK_FNAME(CCTK_ReductionHandle)
     (int *handle, ONE_FORTSTRING_ARG);

void CCTK_FCALL CCTK_FNAME(CCTK_ReductionArrayHandle)
     (int *operation_handle, ONE_FORTSTRING_ARG);

void CCTK_FCALL CCTK_FNAME(CCTK_Reduce)
     (int *fortranreturn,
      const cGH *GH,
      const int *proc,
      const int *operation_handle,
      const int *num_out_vals,
      const int *type_out_vals,
      void *out_vals,
      const int *num_in_fields,
      ... );

void CCTK_FCALL CCTK_FNAME(CCTK_ReduceArray)
     (int *fortran_return,
      const cGH *GH,
      const int *proc,
      const int *operation_handle,
      const int *num_out_vals,
      const int *type_out_vals,
      void *out_vals,
      const int *num_dims,
      const int *num_in_arrays,
      const int *type_in_arrays,
      ... );

/* FIXME: OLD INTERFACE */
void CCTK_FCALL CCTK_FNAME(CCTK_ReduceLocalScalar)
     (int *fortran_return,
      const cGH *GH,
      const int *proc,
      const int *operation_handle,
      const void *in_scalar,
      void *out_scalar,
      const int *data_type);

void CCTK_FCALL CCTK_FNAME(CCTK_ReduceLocScalar)
     (int *fortran_return,
      const cGH *GH,
      const int *proc,
      const int *operation_handle,
      const void *in_scalar,
      void *out_scalar,
      const int *data_type);

void CCTK_FCALL CCTK_FNAME(CCTK_ReduceLocalArray1D)
     (int *fortran_return,
      const cGH *GH,
      const int *proc,
      const int *operation_handle,
      const void *in_array1d,
      void *out_array1d,
      const int *num_in_array1d,
      const int *data_type);

void CCTK_FCALL CCTK_FNAME(CCTK_ReduceLocArrayToArray1D)
     (int *fortran_return,
      const cGH *GH,
      const int *proc,
      const int *operation_handle,
      const void *in_array1d,
      void *out_array1d,
      const int *num_in_array1d,
      const int *data_type);
void CCTK_FCALL  CCTK_FNAME(CCTK_ReduceLocArrayToArray2D)
     (int  *fortran_return, const cGH *GH,
      const int  *proc,
      const int  *operation_handle,
      const void *in_array2d,
      void *out_array2d,
      const int  *xsize, const int *ysize,
      const int  *data_type);
void CCTK_FCALL  CCTK_FNAME(CCTK_ReduceLocArrayToArray3D)
     (int  *fortran_return, const cGH *GH,
      const int  *proc,
      const int  *operation_handle,
      const void *in_array3d,
      void *out_array3d,
      const int  *xsize, const int *ysize, const int *zsize,
      const int  *data_type);
/********************************************************************
 ********************    Internal Typedefs   ************************
 ********************************************************************/
/* structure holding the routines for a registered reduction operator */
typedef struct
{
  const char *implementation;
  const char *name;
  cReduceOperator reduce_operator;
} t_reduce_operator;


/********************************************************************
 ********************    Static Variables   *************************
 ********************************************************************/

static cHandledData *ReductionOperators = NULL;
static int num_reductions = 0;
static cHandledData *ReductionArrayOperators = NULL;
static int num_reductions_array = 0;


 /*@@
   @routine    CCTKi_RegisterReductionOperator
   @date       April 28 1999
   @author     Gabrielle Allen
   @desc
   Registers "function" as a reduction operator called "name"
   @enddesc
   @var     function
   @vdesc   Routine containing reduction operator
   @vtype   (void (*))
   @vio
   @endvar
   @var     name
   @vdesc   String containing name of reduction operator
   @vtype   const char *
   @vio     in
   @endvar
@@*/


int CCTKi_RegisterReductionOperator(const char *thorn,
                                    cReduceOperator operator,
                                    const char *name)
{
  int handle;
  t_reduce_operator *reduce_operator;

  /* Check arguments */

  /* Check that the method hasn't already been registered */
  handle = Util_GetHandle(ReductionOperators, name,
                          (void **) &reduce_operator);

  if(handle < 0)
  {

    reduce_operator = malloc (sizeof (t_reduce_operator));

    if (reduce_operator)
    {
      reduce_operator->implementation =
              CCTK_ThornImplementation(thorn);
      reduce_operator->name = name;
      reduce_operator->reduce_operator = operator;
      handle = Util_NewHandle(&ReductionOperators, name, reduce_operator);

      /* Remember how many reduction operators there are */
      num_reductions++;
    }
  }
  else
  {
    /* Reduction operator with this name already exists. */
    CCTK_Warn(1,__LINE__,__FILE__,"Cactus",
              "CCTK_RegisterReductionOperator: Reduction operator "
              "with this name already exists");
    handle = -1;
  }

#ifdef DEBUG_REDUCTION
  CCTK_PRINTSEPARATOR
  printf("In CCTK_RegisterReductionOperator\n");
  printf("---------------------------------\n");
  printf("  Registering %s with handle %d\n",name,handle);
  CCTK_PRINTSEPARATOR
#endif

  return handle;

}


 /*@@
   @routine    CCTK_ReductionHandle
   @date       April 28 1999
   @author     Gabrielle Allen
   @desc
   Returns the handle of a given reduction operator
   @enddesc
   @var     reduction
   @vdesc   String containing name of reduction operator
   @vtype   const char *
   @vio     in
   @endvar
@@*/

int CCTK_ReductionHandle(const char *reduction)
{

  int handle;
  void **data=NULL; /* isn't used here */

  handle = Util_GetHandle(ReductionOperators, reduction, data);

#ifdef DEBUG_REDUCTION
  CCTK_PRINTSEPARATOR
  printf("In CCTK_ReductionHandle\n");
  printf("-----------------------\n");
  printf("  Got handle %d for %s\n",handle,reduction);
  CCTK_PRINTSEPARATOR
#endif

  if (handle < 0)
    CCTK_VWarn(1,__LINE__,__FILE__,"Cactus",
              "CCTK_ReductionHandle: No handle found for reduction operator '%s'",
              reduction);

  return handle;

}

void CCTK_FCALL CCTK_FNAME(CCTK_ReductionHandle)(int *handle, ONE_FORTSTRING_ARG)
{
  ONE_FORTSTRING_CREATE(reduction)
  *handle = CCTK_ReductionHandle(reduction);
  free(reduction);
}


 /*@@
   @routine    CCTK_Reduce
   @date       April 28 1999
   @author     Gabrielle Allen
   @desc
               Generic routine for doing a reduction operation on a set of
               Cactus variables.
   @enddesc
   @var     GH
   @vdesc   pointer to the grid hierarchy
   @vtype   cGH *
   @vio     in
   @endvar
   @var     proc
   @vdesc   processor that receives the result of the reduction operation
            (a negative value means that all processors get the result)
   @vtype   int
   @vio     in
   @endvar
   @var     operation_handle
   @vdesc   the handle specifying the reduction operator
   @vtype   int
   @vio     in
   @endvar
   @var     num_out_vals
   @vdesc   number of elements in the reduction output
   @vtype   int
   @vio     in
   @endvar
   @var     type_out_vals
   @vdesc   datatype of the output values
   @vtype   int
   @vio     in
   @endvar
   @var     out_vals
   @vdesc   pointer to buffer holding the output values
   @vtype   void *
   @vio     in
   @endvar
   @var     num_in_fields
   @vdesc   number of input fields passed in the variable argument list
   @vtype   int
   @vio     in
   @endvar
   @var     <...>
   @vdesc   list of variables indices of input fields
   @vtype   int
   @vio     in
   @endvar
@@*/

int CCTK_Reduce(const cGH *GH,
                int proc,
                int operation_handle,
                int num_out_vals,
                int type_out_vals,
                void *out_vals,
                int num_in_fields,
                ... )
{
  va_list indices;
  int i;
  int retval;
  int *in_fields;
  t_reduce_operator *operator;

  /* Get the pointer to the reduction operator */
  if (operation_handle < 0)
  {
    CCTK_Warn(3,__LINE__,__FILE__,"Cactus",
              "CCTK_Reduce: Invalid handle passed to CCTK_Reduce");
    retval = -1;
  }
  else
  {
    operator = Util_GetHandledData(ReductionOperators,operation_handle);

    if (!operator)
    {
      CCTK_Warn(3,__LINE__,__FILE__,"Cactus",
                "CCTK_Reduce: Reduction operation is not registered"
                "and cannot be called");
      retval = -1;
    }
    else
    {
      /* Fill in the array of variable indices from the variable
         argument list */
      in_fields = malloc(num_in_fields*sizeof(int));
      va_start(indices, num_in_fields);
      for (i=0; i<num_in_fields; i++)
      {
        in_fields[i] = va_arg(indices,int);
      }
      va_end(indices);

      retval = operator->reduce_operator (GH, proc, num_out_vals,
                                        type_out_vals, out_vals,
                                        num_in_fields, in_fields);

      free(in_fields);
    }
  }
  return retval;
}

void CCTK_FCALL CCTK_FNAME(CCTK_Reduce)
     (int *fortranreturn,
      const cGH *GH,
      const int *proc,
      const int *operation_handle,
      const int *num_out_vals,
      const int *type_out_vals,
      void *out_vals,
      const int *num_in_fields,
      ... )
{
  va_list indices;
  int retval;
  int i;
  int *in_fields;
  t_reduce_operator *operator;

  /* initialize return code to indicate an error */
  *fortranreturn = -1;

  if (*operation_handle < 0)
  {
    CCTK_Warn(3,__LINE__,__FILE__,"Cactus",
              "CCTK_Reduce: Invalid handle passed to CCTK_Reduce");
    retval = -1;
  }
  else
  {
    /* Get the pointer to the reduction operator */
    operator = Util_GetHandledData(ReductionOperators,*operation_handle);

    if (!operator)
    {
      CCTK_Warn(3,__LINE__,__FILE__,"Cactus",
                "CCTK_Reduce: Reduction operation is not registered"
                " and cannot be called");
      retval = -1;
    }
    else
    {

      /* Fill in the array of variable indices from the variable
         argument list */
      in_fields = malloc (*num_in_fields * sizeof (int));
      va_start(indices, num_in_fields);
      for (i=0; i<*num_in_fields; i++)
      {
        in_fields[i] = *va_arg(indices,int *);
      }
      va_end(indices);

      retval = operator->reduce_operator (GH, *proc, *num_out_vals,
                                          *type_out_vals, out_vals,
                                          *num_in_fields,in_fields);

      free(in_fields);
    }
  }
  *fortranreturn = retval;
}


 /*@@
   @routine CCTK_RegisterReductionArrayOperator
   @date    Aug 19 1999
   @author  Thomas Radke
   @desc
            Registers "function" as a array reduction operator called "name"
   @enddesc
   @var     function
   @vdesc   Routine containing reduction operator
   @vtype   (int (*))
   @vio
   @endvar
   @var     name
   @vdesc   String containing name of reduction operator
   @vtype   const char *
   @vio     in
   @endvar
@@*/

int CCTK_RegisterReductionArrayOperator
         (int (*function)(REDUCTION_ARRAY_OPERATOR_REGISTER_ARGLIST),
         const char *name)
{
  int handle;

  /* Check that the method hasn't already been registered */
  handle = Util_GetHandle(ReductionArrayOperators, name, NULL);

  if(handle < 0)
  {
    /* Get a handle for it. */
    handle = Util_NewHandle(&ReductionArrayOperators, name, (int *)function);

    /* Remember how many reduction operators there are */
    num_reductions_array++;
   }
  else
  {
    /* Reduction operator with this name already exists. */
    CCTK_Warn(1,__LINE__,__FILE__,"Cactus",
              "CCTK_RegisterReductionArrayOperator: "
              "Array reduction operator with this name already exists");
    handle = -1;
  }

#ifdef DEBUG_REDUCTION
  CCTK_PRINTSEPARATOR
  printf("In CCTK_RegisterReductionArrayOperator\n");
  printf("---------------------------------\n");
  printf("  Registering %s with handle %d\n",name,handle);
  CCTK_PRINTSEPARATOR
#endif

  return handle;

}


 /*@@
   @routine CCTK_ReductionArrayHandle
   @date    Aug 19 1999
   @author  Thomas Radke
   @desc
            Returns the handle of a given array reduction operator
   @enddesc
   @var     reduction
   @vdesc   String containing name of array reduction operator
   @vtype   const char *
   @vio     in
   @endvar
@@*/

int CCTK_ReductionArrayHandle(const char *reduction)
{

  int handle;
  void **data=NULL; /* isn't used here */

  handle = Util_GetHandle(ReductionArrayOperators, reduction, data);

#ifdef DEBUG_REDUCTION
  CCTK_PRINTSEPARATOR
  printf("In CCTK_ReductionArrayHandle\n");
  printf("-----------------------\n");
  printf("  Got handle %d for %s\n",handle,reduction);
  CCTK_PRINTSEPARATOR
#endif

  if (handle < 0)
  {
    CCTK_VWarn (1, __LINE__, __FILE__, "Cactus",
                "CCTK_ReductionArrayHandle: "
                "No handle found for array reduction operator '%s'",
                reduction);
  }

  return handle;
}

void CCTK_FCALL CCTK_FNAME(CCTK_ReductionArrayHandle)
     (int *operation_handle, ONE_FORTSTRING_ARG)
{
  ONE_FORTSTRING_CREATE(reduction)
  *operation_handle = CCTK_ReductionArrayHandle(reduction);
  free(reduction);
}


 /*@@
   @routine    CCTK_ReduceArray
   @date       Aug 19 1999
   @author     Thomas Radke
   @desc
               Generic routine for doing a reduction operation on a set of
               arrays.
   @enddesc
   @var     GH
   @vdesc   pointer to the grid hierarchy
   @vtype   const cGH *
   @vio     in
   @endvar
   @var     proc
   @vdesc   processor that receives the result of the reduction operation
            (a negative value means that all processors get the result)
   @vtype   int
   @vio     in
   @endvar
   @var     operation_handle
   @vdesc   the handle specifying the reduction operator
   @vtype   int
   @vio     in
   @endvar
   @var     num_out_vals
   @vdesc   number of elements in the reduction output
   @vtype   int
   @vio     in
   @endvar
   @var     type_out_vals
   @vdesc   datatype of the reduction output
   @vtype   int
   @vio     in
   @endvar
   @var     out_vals
   @vdesc   pointer to buffer holding the output values
   @vtype   void *
   @vio     in
   @endvar
   @var     num_dims
   @vdesc   number of dimensions of input arrays
   @vtype   int
   @vio     in
   @endvar
   @var     num_in_fields
   @vdesc   number of input fields passed in the variable argument list
   @vtype   int
   @vio     in
   @endvar
   @var     type_in_arrays
   @vdesc   datatype of the input arrays
   @vtype   int
   @vio     in
   @endvar
   @var     <...>
   @vdesc   list of dimensions of input arrays
   @vtype   int
   @vio     in
   @endvar
   @var     <...>
   @vdesc   list of pointers to input arrays
   @vtype   void *
   @vio     in
   @endvar
@@*/

int CCTK_ReduceArray(const cGH *GH,
                     int proc,
                     int operation_handle,
                     int num_out_vals,
                     int type_out_vals,
                     void *out_vals,
                     int num_dims,
                     int num_in_arrays,
                     int type_in_arrays,
                     ... )
{

  va_list indices;
  int i;
  int *dims;
  const void **in_arrays;
  int (*function)(REDUCTION_ARRAY_OPERATOR_REGISTER_ARGLIST)=NULL;


  /* Get the pointer to the reduction operator */
  if (operation_handle < 0)
  {
    CCTK_Warn(3,__LINE__,__FILE__,"Cactus",
              "CCTK_ReduceArray: Invalid handle passed to CCTK_ReduceArray");
    return (-1);
  }

  function = Util_GetHandledData(ReductionArrayOperators,operation_handle);

  if (! function)
  {
    CCTK_Warn(3,__LINE__,__FILE__,"Cactus",
              "CCTK_ReduceArray: Array reduction operation is not registered "
                 "and cannot be called");
    return (-1);
  }

  /* allocate memory for dims and input array pointers */
  dims = malloc (num_dims * sizeof (int));
  in_arrays = malloc (num_in_arrays * sizeof (void *));

  /* Fill in the arrays of dims and input array pointers
     from the variable argument list */

  va_start(indices, type_in_arrays);

  for (i = 0; i < num_dims; i++)
  {
    dims[i] = va_arg (indices, int);
  }
  for (i = 0; i < num_in_arrays; i++)
  {
    in_arrays[i] = va_arg (indices, void *);
  }

  va_end(indices);

  function (GH, proc, num_dims, dims,
            num_in_arrays, in_arrays, type_in_arrays,
            num_out_vals, out_vals, type_out_vals);

  free (in_arrays);
  free (dims);

  return (0);
}

void CCTK_FCALL CCTK_FNAME(CCTK_ReduceArray)
     (int *fortran_return,
      const cGH *GH,
      const int *proc,
      const int *operation_handle,
      const int *num_out_vals,
      const int *type_out_vals,
      void *out_vals,
      const int *num_dims,
      const int *num_in_arrays,
      const int *type_in_arrays,
      ... )
{

  va_list varargs;
  int i;
  int *dims;
  const void **in_arrays;
  int (*function)(REDUCTION_ARRAY_OPERATOR_REGISTER_ARGLIST);


  /* initialize return code to indicate an error */
  *fortran_return = -1;

  /* Get the pointer to the reduction operator */
  if (*operation_handle < 0)
  {
    CCTK_Warn (3,__LINE__,__FILE__,"Cactus",
               "CCTK_ReduceArray: Invalid handle passed to CCTK_ReduceArray");
    return;
  }

  function = Util_GetHandledData (ReductionArrayOperators, *operation_handle);

  if (! function)
  {
    CCTK_Warn (3,__LINE__,__FILE__,"Cactus",
               "CCTK_ReduceArray: Array reduction operation is not registered "
               "and cannot be called");
    return;
  }

  /* allocate memory for dims and input array pointers */
  dims = malloc (*num_dims * sizeof (int));
  in_arrays = malloc (*num_in_arrays * sizeof (void *));

  /* Fill in the arrays of dims and input array pointers
     from the variable argument list */

  va_start (varargs, type_in_arrays);

  for (i = 0; i < *num_dims; i++)
  {
    dims[i] = *va_arg (varargs, int *);
  }
  for (i = 0; i < *num_in_arrays; i++)
  {
    in_arrays[i] = va_arg (varargs, void *);
  }

  va_end (varargs);

  function (GH, *proc, *num_dims, dims,
            *num_in_arrays, in_arrays, *type_in_arrays,
            *num_out_vals, out_vals, *type_out_vals);

  free (in_arrays);
  free (dims);

  *fortran_return = 0;
}


 /*@@
   @routine    CCTK_ReduceLocalScalar
   @date       Aug 19 1999
   @author     Thomas Radke
   @desc    Wrapper function for reduction of a single scalar
   @enddesc
   @var     GH
   @vdesc   pointer to the grid hierarchy
   @vtype   const cGH *
   @vio     in
   @endvar
   @var     proc
   @vdesc   processor that receives the result of the reduction operation
            (a negative value means that all processors get the result)
   @vtype   int
   @vio     in
   @endvar
   @var     operation_handle
   @vdesc   the handle specifying the reduction operator
   @vtype   int
   @vio     in
   @endvar
   @var     in_scalar
   @vdesc   pointer to input scalar
   @vtype   void *
   @vio     in
   @endvar
   @var     out_scalar
   @vdesc   pointer to output scalar
   @vtype   void *
   @vio     in
   @endvar
   @var     data_type
   @vdesc   datatype for both input and output scalar
   @vtype   int
   @vio     in
   @endvar
@@*/

/*** FIXME: OLD INTERFACE gerd ***/
int CCTK_ReduceLocalScalar (const cGH *GH, int proc, int operation_handle,
                            const void *in_scalar, void *out_scalar, int data_type)
{
  return (CCTK_ReduceArray (GH, proc, operation_handle,
                            1, data_type, out_scalar,
                            1, 1, data_type, 1, in_scalar));
}

/*** FIXME: OLD INTERFACE gerd ***/
void CCTK_FCALL CCTK_FNAME(CCTK_ReduceLocalScalar)
     (int *fortran_return,
      const cGH *GH,
      const int *proc,
      const int *operation_handle,
      const void *in_scalar,
      void *out_scalar,
      const int *data_type)
{
  *fortran_return = CCTK_ReduceArray (GH, *proc, *operation_handle,
                                      1, *data_type, out_scalar,
                                      1, 1, *data_type, 1, in_scalar);
}


int CCTK_ReduceLocScalar (const cGH *GH, int proc, int operation_handle,
                            const void *in_scalar, void *out_scalar, int data_type)
{
  return (CCTK_ReduceArray (GH, proc, operation_handle,
                            1, data_type, out_scalar,
                            1, 1, data_type, 1, in_scalar));
}

void CCTK_FCALL CCTK_FNAME(CCTK_ReduceLocScalar)
     (int *fortran_return,
      const cGH *GH,
      const int *proc,
      const int *operation_handle,
      const void *in_scalar,
      void *out_scalar,
      const int *data_type)
{
  *fortran_return = CCTK_ReduceArray (GH, *proc, *operation_handle,
                                      1, *data_type, out_scalar,
                                      1, 1, *data_type, 1, in_scalar);
}


 /*@@
   @routine    CCTK_ReduceLocArrayToArray1D
   @date       Thu Oct 14 12:10:01 1999
   @author     Gerd Lanfermann
   @desc
        Interface to the migthy CCTK_Reduce for
        reduction of local 1D arrays to local arrays
        (element by element).
   @enddesc
@@*/

/*** FIXME: OLD INTERFACE gerd ***/
int CCTK_ReduceLocalArray1D (const cGH *GH, int proc, int operation_handle,
                             const void *in_array1d, void *out_array1d,
                             int num_in_array1d,
                             int data_type)
{
  return (CCTK_ReduceArray (GH, proc, operation_handle,
                            num_in_array1d, data_type, out_array1d,
                            1, 1, data_type, num_in_array1d, in_array1d));
}

int CCTK_ReduceLocArrayToArray1D(const cGH *GH, int proc, int operation_handle,
                                 const void *in_array1d, void *out_array1d,
                                 int num_in_array1d,
                                 int data_type)
{
  return (CCTK_ReduceArray (GH, proc, operation_handle,
                            num_in_array1d, data_type, out_array1d,
                            1, 1, data_type, num_in_array1d, in_array1d));
}

/*** FIXME: OLD INTERFACE gerd ***/
void CCTK_FCALL CCTK_FNAME(CCTK_ReduceLocalArray1D)
     (int *fortran_return,
      const cGH *GH,
      const int *proc,
      const int *operation_handle,
      const void *in_array1d,
      void *out_array1d,
      const int *num_in_array1d,
      const int *data_type)
{
  *fortran_return = CCTK_ReduceArray (GH, *proc, *operation_handle,
                                      *num_in_array1d, *data_type, out_array1d,
                                      1, 1, *data_type, *num_in_array1d,
                                      in_array1d);
}


/*@@
   @routine    CCTK_ReduceLocArrayToArray1D
   @date       Sat Nov 27 22:52:10 1999
   @author     Gerd Lanfermann
   @desc
       Interface for the reduction of local 1d arrays
       to the mighty CCCTK_reduce interface.
   @enddesc
@@*/


void CCTK_FCALL CCTK_FNAME(CCTK_ReduceLocArrayToArray1D)
     (int *fortran_return,
      const cGH *GH,
      const int *proc,
      const int *operation_handle,
      const void *in_array1d,
      void *out_array1d,
      const int *num_in_array1d,
      const int *data_type)
{
  *fortran_return = CCTK_ReduceArray (GH, *proc, *operation_handle,
                                      *num_in_array1d, *data_type, out_array1d,
                                      1, 1, *data_type, *num_in_array1d,
                                      in_array1d);
}

/*@@
   @routine    CCTK_ReduceLocArrayToArray2D
   @date       Sat Nov 27 22:52:10 1999
   @author     Gerd Lanfermann
   @desc
       Interface for the reduction of local 2d arrays
       to the mighty CCCTK_reduce interface.
   @enddesc
@@*/



int CCTK_ReduceLocArrayToArray2D(const cGH *GH, int proc, int operation_handle,
                                 const void *in_array2d, void *out_array2d,
                                 int xsize, int ysize,
                                 int data_type)
{
  int  lin_size= xsize*ysize;
  return (CCTK_ReduceArray (GH, proc, operation_handle,
                            lin_size,
                            data_type, out_array2d,
                            2, 1, data_type,
                            xsize,ysize, in_array2d));
}

void CCTK_FCALL  CCTK_FNAME(CCTK_ReduceLocArrayToArray2D)
     (int  *fortran_return, const cGH *GH,
      const int  *proc,
      const int  *operation_handle,
      const void *in_array2d,
      void *out_array2d,
      const int  *xsize, const int *ysize,
      const int  *data_type)
{
  int lin_size = (*xsize)*(*ysize);
  *fortran_return =  CCTK_ReduceArray (GH, *proc, *operation_handle,
                                      lin_size,
                                      *data_type, out_array2d,
                                      2, 1, *data_type,
                                      *xsize, *ysize,
                                      in_array2d);
}

/*@@
   @routine    CCTK_ReduceLocArrayToArray1D
   @date       Sat Nov 27 22:52:10 1999
   @author     Gerd Lanfermann
   @desc
       Interface for the reduction of local 3d arrays
       to 3d arrays.
   @enddesc
@@*/

int CCTK_ReduceLocArrayToArray3D(const cGH *GH, int proc, int operation_handle,
                                 const void *in_array3d, void *out_array3d,
                                 int xsize, int  ysize, int zsize,
                                 int data_type)
{

  int lin_size =  xsize*ysize*zsize;
  return (CCTK_ReduceArray (GH, proc, operation_handle,
                            lin_size,
                            data_type, out_array3d,
                            3, 1, data_type,
                            xsize,ysize,zsize,
                            in_array3d));
}

void CCTK_FCALL  CCTK_FNAME(CCTK_ReduceLocArrayToArray3D)
     (int  *fortran_return, const cGH *GH,
      const int  *proc,
      const int  *operation_handle,
      const void *in_array3d,
      void *out_array3d,
      const int  *xsize, const int *ysize, const int *zsize,
      const int  *data_type)
{
  int lin_size =  (*xsize)*(*ysize)*(*zsize);
  *fortran_return =  CCTK_ReduceArray (GH, *proc, *operation_handle,
                                       lin_size,
                                       *data_type, out_array3d,
                                       3, 1, *data_type,
                                       *xsize,*ysize,*zsize,
                                       in_array3d);
}

 /*@@
   @routine    CCTK_NumReduceOperators
   @date       Mon Oct 22 2001
   @author     Gabrielle Allen
   @desc
               The number of reduction operators registered
   @enddesc
   @returntype int
   @returndesc
               number of reduction operators
   @endreturndesc
@@*/

int CCTK_NumReduceOperators()
{
  return num_reductions;
}

 /*@@
   @routine    CCTK_ReduceOperatorImplementation
   @date       Mon Oct 22 2001
   @author     Gabrielle Allen
   @desc
   Provide the implementation which provides an reduction operator
   @enddesc
   @returntype int
   @returndesc
               Implementation which supplied the interpolation operator
   @endreturndesc
@@*/

const char *CCTK_ReduceOperatorImplementation(int handle)
{
  t_reduce_operator *operator;
  const char *imp=NULL;

  operator = Util_GetHandledData (ReductionOperators, handle);

  if (operator)
  {
    imp = operator->implementation;
  }

  return imp;
}



 /*@@
   @routine    CCTK_ReduceOperator
   @date       December 27 2001
   @author     Gabrielle Allen
   @desc
               Returns the name of a reduction operator
   @enddesc
   @var        handle
   @vdesc      Handle for reduction operator
   @vtype      int
   @vio        in
   @endvar

   @returntype const char *
   @returndesc
   The name of the reduction operator, or NULL if the handle
   is invalid
   @endreturndesc
@@*/
const char *CCTK_ReduceOperator (int handle)
{
  const char *name=NULL;
  t_reduce_operator *operator;

  if (handle < 0)
  {
    CCTK_VWarn (6, __LINE__, __FILE__, "Cactus",
                "CCTK_ReduceOperator: Handle %d invalid", handle);
  }
  else
  {
    operator = Util_GetHandledData (ReductionOperators, handle);
    if (operator)
    {
      name = operator->name;
    }
    else
    {
      CCTK_VWarn (6, __LINE__, __FILE__, "Cactus",
                  "CCTK_ReduceOperator: Handle %d invalid", handle);
    }
  }

  return name;
}