aboutsummaryrefslogtreecommitdiff
path: root/src/Comm.c
blob: 792882e3a2731dbf0577e1db4041dc297578e30b (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
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
 /*@@
   @file      Comm.c
   @date      Thu Feb  4 11:34:29 1999
   @author    Tom Goodale
   @desc
              PUGH communication functions
   @enddesc
   @version   $Id$
 @@*/

/*#define DEBUG_PUGH 1*/

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

#include "cctk.h"

#include "pugh.h"
#include "pughi.h"
#include "pugh_Comm.h"

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

CCTK_FILEVERSION(CactusPUGH_PUGH_Comm_c);


/* local function prototypes */
static int PUGH_EnableGArrayGroupComm(pGH *pughGH, int first_var, int commflag);
static int PUGH_EnableComm(pGH *pughGH, pComm *comm, int commflag);
static int PUGH_DisableComm(pGA *GA, pComm *comm);
static int PUGH_SyncGArrayGroup(pGH *pughGH, int first_var);
static int PUGH_Sync(pGH *pughGH, pComm *comm);
static int PUGH_SyncSingleProc(pGH *pughGH, pComm *comm);


/*@@
   @routine    PUGH_SyncGroup
   @author     Thomas Radke
   @date       30 Mar 1999
   @desc
               Synchronizes all variables in the group indicated by groupname.
   @enddesc
   @calls      PUGH_SyncGArrayGroup

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        groupname
   @vdesc      name of the group to be synchronized
   @vtype      const char *
   @vio        in
   @endvar

   @returntype int
   @returndesc
                  return code of @seeroutine PUGH_SyncGArrayGroup
                  for grid array groups<BR>
               -1 for an unknown group or group type
   @endreturndesc
@@*/
int PUGH_SyncGroup(const cGH *GH, const char *groupname)
{
  cGroup pgroup;
  int group, retval;


#ifdef DEBUG_PUGH
  printf (" PUGH_SyncGroup: request for group '%s'\n", groupname);
  fflush (stdout);
#endif

  /* get the group info from its index */
  group = CCTK_GroupIndex(groupname);
  if (group < 0)
  {
    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                "PUGH_SyncGroup: Unknown group: %s", groupname);
    retval = -1;
  }
  else if (CCTK_NumVarsInGroupI(group) == 0)
  {
    retval = 0;  /* do nothing */
  }
  else
  {
    CCTK_GroupData(group, &pgroup);

    if (pgroup.numvars == 0)
    {
      retval = 0;
    }
    else
    {
      retval = PUGH_SyncGArrayGroup(PUGH_pGH(GH), CCTK_FirstVarIndexI(group));
    }
  }

  return (retval);
}


/*@@
   @routine    PUGH_SyncGroupsByDirI
   @author     Thomas Radke
   @date       29 June 2006
   @desc
               Synchronizes a list of groups in the given directions.
   @enddesc
   @calls      PUGH_SyncGArrayGroup

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        num_groups
   @vdesc      number of groups to be synchronised
   @vtype      int
   @vio        in
   @endvar
   @var        groups
   @vdesc      list of indices of groups to be synchronised
   @vtype      const int *
   @vio        in
   @endvar
   @var        directions
   @vdesc      list of directions to be synchronised
   @vtype      const int *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               returns the number of groups synchronised successfully
   @endreturndesc
@@*/
int PUGH_SyncGroupsByDirI(const cGH *GH,
                          int num_groups,
                          const int *groups,
                          const int *directions)
{
  int group, first_var, retval = 0;
  pGH *pughGH = PUGH_pGH (GH); 


  /* individual directions aren't supported (yet?) */
  if (directions != NULL)
  {
    CCTK_WARN (0, "PUGH doesn't support synchronisation of individual "
                  "directions");
  }

  /* synchronise all groups one by one */
  for (group = 0; group < num_groups; group++)
  {
#ifdef DEBUG_PUGH
    char *groupname = CCTK_GroupName (groups[group]);
    printf (" PUGH_SyncGroupsByDirI: request for group '%s'\n", groupname);
    fflush (stdout);
    free (groupname);
#endif

    first_var = CCTK_FirstVarIndexI (groups[group]);
    if (first_var >= 0 && PUGH_SyncGArrayGroup (pughGH, first_var) == 0)
    {
      retval++;
    }
  }

  return (retval);
}


/*@@
   @routine PUGH_EnableGroupComm
   @author     Thomas Radke
   @date       30 Mar 1999
   @desc
               Enables communication for all variables in the group
               indicated by groupname.
   @enddesc
   @calls      PUGH_EnableGArrayGroupComm

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        groupname
   @vdesc      name of the group to be synchronized
   @vtype      const char *
   @vio        in
   @endvar

   @returntype int
   @returndesc
                  return code of @seeroutine PUGH_EnableGArrayGroupComm
                  for grid array groups<BR>
               -1 for an unknown group type
   @endreturndesc
 @@*/
int PUGH_EnableGroupComm(const cGH *GH, const char *groupname)
{
  int group, retval;
  cGroup pgroup;


#ifdef DEBUG_PUGH
  printf(" PUGH_EnableGroupComm: request for group '%s'\n", groupname);
  fflush(stdout);
#endif

  /* get the group info from its index */
  group = CCTK_GroupIndex(groupname);
  CCTK_GroupData(group, &pgroup);

  retval = PUGH_EnableGArrayGroupComm(PUGH_pGH(GH),
                                      CCTK_FirstVarIndexI(group),
                                      PUGH_ALLCOMM);

  return (retval);
}


/*@@
   @routine    PUGH_DisableGroupComm
   @author     Thomas Radke
   @date       30 Mar 1999
   @desc
               Disables communication for all variables in the group
               indicated by groupname.
   @enddesc
   @calls      PUGH_DisableGArrayGroupComm

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        groupname
   @vdesc      name of the group to be synchronized
   @vtype      const char *
   @vio        in
   @endvar

   @returntype int
   @returndesc
                  return code of @seeroutine PUGH_DisableGArrayGroupComm
                  for grid array groups<BR>
               -1 for an unknown group type
   @endreturndesc
 @@*/
int PUGH_DisableGroupComm(const cGH *GH, const char *groupname)
{
  int group, first_var, retval;
  cGroup pgroup;
  pGA *GA;
  pGH *pughGH;


#ifdef DEBUG_PUGH
  printf(" PUGH_DisableGroupComm: request for group '%s'\n", groupname);
  fflush(stdout);
#endif

  /* get the group info from its index */
  group = CCTK_GroupIndex(groupname);
  CCTK_GroupData(group, &pgroup);

  first_var = CCTK_FirstVarIndexI(group);
  if (first_var < 0)
  {
    CCTK_WARN (0, "Illegal group index -- group has no variables");
  }
  pughGH = PUGH_pGH(GH);
  GA = (pGA *) pughGH->variables[first_var][0];

  /* FIXME:  workaround.  This one is really bad ! */
  retval = PUGH_DisableGArrayGroupComm(pughGH, first_var, GA->groupcomm);

  return (retval);
}


/*@@
   @routine    PUGH_EnableGArrayComm
   @date       Mon Jun 05 2000
   @author     Thomas Radke
   @desc
               Enables communication for a single array.
   @enddesc

   @var        GA
   @vdesc      Pointer to grid array structure
   @vtype      pGA *
   @vio        in
   @endvar
   @var        commflag
   @vdesc      flag indicating in which directions to communicate
   @vtype      int
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine PUGH_EnableComm
   @endreturndesc
 @@*/
int PUGH_EnableGArrayComm(pGA *GA, int commflag)
{

#ifdef DEBUG_PUGH
  printf(" PUGH_EnableGArrayComm: request for var '%s' commflag %d\n",
          GA->name, commflag);
  fflush(stdout);
#endif

  return (PUGH_EnableComm((pGH *) GA->parent, GA->comm, commflag));
}


/*@@
   @routine    PUGH_DisableGArrayComm
   @date       Mon Jun 05 2000
   @author     Thomas Radke
   @desc
               Disables communication for a single array.
   @enddesc

   @var        GA
   @vdesc      Pointer to grid array structure
   @vtype      pGA *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine PUGH_DisableComm
   @endreturndesc
 @@*/
int PUGH_DisableGArrayComm(pGA *GA)
{
#ifdef DEBUG_PUGH
  printf(" PUGH_DisableGArrayComm: request for var '%s'\n", GA->name);
  fflush(stdout);
#endif

  return (PUGH_DisableComm(GA, GA->comm));
}


/*@@
   @routine    PUGH_SyncGArrayGroup
   @date       Mon Jun 05 2000
   @author     Thomas Radke
   @desc
               Synchronizes a group of arrays
               given the first variable within this group.
   @enddesc

   @var        pughGH
   @vdesc      Pointer to PUGH grid hierarchy extension
   @vtype      pGH *
   @vio        in
   @endvar
   @var        first_var
   @vdesc      grid index of first variable in the group
   @vtype      int
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine PUGH_Sync
   @endreturndesc
 @@*/
static int PUGH_SyncGArrayGroup(pGH *pughGH, int first_var)
{
  int retval;
  pGA *firstGA;


  if (first_var<0 || first_var>=CCTK_NumVars())
  {
    CCTK_WARN (0, "internal error");
  }

  firstGA = (pGA *) pughGH->variables [first_var][0];

#ifdef DEBUG_PUGH
  printf(" PUGH_SyncGArrayGroup: request for group with first var '%s'\n",
         firstGA->name);
  fflush(stdout);
#endif

  /* Automatically enable communication for this group. */
  PUGH_EnableGArrayGroupComm(pughGH, first_var, PUGH_ALLCOMM);
  retval = PUGH_Sync(pughGH, firstGA->groupcomm);

  return (retval);
}


/*@@
   @routine    PUGH_SyncGArray
   @date       Mon Jun 05 2000
   @author     Thomas Radke
   @desc
               Synchronizes a single array variable given the GA structure
               of this variable.
   @enddesc

   @var        GA
   @vdesc      Pointer to grid array structure
   @vtype      pGA *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine PUGH_Sync
   @endreturndesc
 @@*/
int PUGH_SyncGArray(pGA *GA)
{

#ifdef DEBUG_PUGH
  printf(" PUGH_SyncGArray: request for var '%s'\n", GA->name);
  fflush(stdout);
#endif

  return (PUGH_Sync((pGH *) GA->parent, GA->comm));
}


/*@@
   @routine    PUGH_DisableGArrayGroupComm
   @date       Mon Jun 05 2000
   @author     Thomas Radke
   @desc
               Disables communication for a group of arrays
               given the first variable within this group.
   @enddesc

   @var        pughGH
   @vdesc      Pointer to PUGH grid hierarchy extension
   @vtype      pGH *
   @vio        in
   @endvar
   @var        first_var
   @vdesc      grid index of first variable in the group
   @vtype      int
   @vio        in
   @endvar
   @var        groupcomm
   @vdesc      pointer to comm structure for this group
   @vtype      pComm *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine PUGH_DisableComm
   @endreturndesc
@@*/
int PUGH_DisableGArrayGroupComm(pGH *pughGH, int first_var, pComm *groupcomm)
{
  pGA *GA;                 /* first variable in group */


  if (first_var < 0)
  {
    CCTK_WARN (0, "Illegal variable index");
  }

  GA = (pGA *) pughGH->variables[first_var][0];

#ifdef DEBUG_PUGH
  printf(" PUGH_DisableGArrayGroupComm: request for group "
         "with first var '%s'\n", GA->name);
  fflush(stdout);
#endif

  return (PUGH_DisableComm(GA, groupcomm));
}


/*@@
   @routine    PUGH_Barrier
   @date       Mon Jun 05 2000
   @author     Thomas Radke
   @desc
               Synchronize all processors at a given point of execution.
   @enddesc

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               0 for success
   @endreturndesc
 @@*/
int PUGH_Barrier(const cGH *GH)
{
#ifdef CCTK_MPI
  CACTUS_MPI_ERROR (MPI_Barrier (GH ? PUGH_pGH (GH)->PUGH_COMM_WORLD :
                                      MPI_COMM_WORLD));
#else
  GH = GH;
#endif

  return (0);
}


/*****************************************************************************/
/* local functions                                                           */
/*****************************************************************************/

/*@@
   @routine    PUGH_EnableGArrayGroupComm
   @date       Mon Jun 05 2000
   @author     Thomas Radke
   @desc
               Enables communication for a group of arrays
               given the first variable within this group.
   @enddesc

   @var        pughGH
   @vdesc      Pointer to PUGH grid hierarchy extension
   @vtype      pGH *
   @vio        in
   @endvar
   @var        first_var
   @vdesc      grid index of first variable in the group
   @vtype      int
   @vio        in
   @endvar
   @var        commflag
   @vdesc      flag indicating in which directions to communicate
   @vtype      int
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine PUGH_EnableComm
   @endreturndesc
 @@*/
static int PUGH_EnableGArrayGroupComm(pGH *pughGH, int first_var, int commflag)
{
  pGA *GA;                 /* first variable in group */


  if (first_var < 0)
  {
    CCTK_WARN (0, "Illegal variable index");
  }

  GA = pughGH->variables [first_var][0];

#ifdef DEBUG_PUGH
  printf(" PUGH_EnableGArrayGroupComm: request for group "
         "with first var '%s', commflag %d\n", GA->name, commflag);
  fflush(stdout);
#endif

  return (PUGH_EnableComm(pughGH, GA->groupcomm, commflag));
}


/*@@
   @routine    PUGH_EnableComm
   @date       Sun Jan 23 12:46:23 2000
   @author     Gabrielle Allen
   @desc
               This sets the docomm[2*dim] array of the GA based
               on the setting of the comm flag and allocates the comm buffers.
   @enddesc
   @history
               @date Mon Jun 05 2000 @author Thomas Radke
               Moved buffer allocation from PUGH_EnableGArrayDataStorage
   @endhistory

   @var        pughGH
   @vdesc      Pointer to PUGH grid hierarchy extension
   @vtype      pGH *
   @vio        in
   @endvar
   @var        comm
   @vdesc      Pointer to comm structure
   @vtype      pComm *
   @vio        in
   @endvar
   @var        commflag
   @vdesc      flag indicating in which directions to communicate
   @vtype      int
   @vio        in
   @endvar

   @returntype int
   @returndesc
               1 if communication was enabled,<BR>
               or negative otherwise
   @endreturndesc
@@*/
static int PUGH_EnableComm(pGH *pughGH, pComm *comm, int commflag)
{
  int retval;              /* return value */
#ifdef CCTK_MPI
  pGA *GA;                 /* GA structure the pComm belongs to */
  int i, idir;             /* loopers */
  int dir;                 /* direction */
  int sz;                  /* buffer size */
#endif


  retval = 1;

#ifdef CCTK_MPI
  /* check whether this comm buffer is already set up */
  if (comm->commflag == commflag)
  {
    return (retval);
  }

  /* get the GA structure the comm structure belongs to
     For a group comm structure this GA is the first one within the group. */
  GA = (pGA *) pughGH->variables [comm->first_var][comm->sync_timelevel];

  if (comm->commflag == PUGH_NOCOMM)
  {

#ifdef DEBUG_PUGH
    printf (" PUGH_EnableComm: allocating comm buffer for %d vars starting "
            "with %d '%s'\n", comm->n_vars, comm->first_var, GA->name);
    fflush (stdout);
#endif

    /* allocate memory for communication buffers: 2 faces per direction */
    for (i = 0; i < 2 * GA->extras->dim; i++)
    {
      comm->buffer_sz[i]   = 0;
      comm->send_buffer[i] = NULL;
      comm->recv_buffer[i] = NULL;

      if (GA->connectivity->neighbours[pughGH->myproc][i] >= 0)
      {
        dir = i/2;

        /* no ghostzones -> no comm buffers */
        sz = comm->n_vars * GA->extras->nghostzones[dir];
        if (sz > 0)
        {
          for (idir = 0; idir < GA->extras->dim; idir++)
          {
            if (idir != dir)
            {
              sz *= GA->extras->lnsize[idir];
            }
          }
          comm->buffer_sz[i]   = sz;
          comm->send_buffer[i] = malloc(sz * GA->varsize);
          comm->recv_buffer[i] = malloc(sz * GA->varsize);

          if (! (comm->send_buffer[i] && comm->recv_buffer[i]))
          {
            for (; i >=0 ; i--)
            {
              free(comm->send_buffer[i]);
              free(comm->recv_buffer[i]);
              comm->buffer_sz[i] = 0;
            }

            CCTK_WARN(1, "Out of memory !");
            retval = -1;
            break;
          }
        }
      }
    }
  }

  /* set up comm flags for each face */
  if (retval >= 0)
  {
    /* Copy commflag */
    comm->commflag = commflag;

    /* First set all communcation off */
    for (idir = 0; idir < 2 * GA->extras->dim; idir++)
    {
      comm->docomm[idir] = 0;
    }

    if (commflag == PUGH_ALLCOMM)
    {
      for (idir = 0; idir < 2 * GA->extras->dim; idir++)
      {
        comm->docomm[idir] = comm->buffer_sz[idir] > 0;
      }
    }
    else if (commflag == PUGH_PLUSFACESCOMM)
    {
      for (idir = 0; idir < GA->extras->dim; idir++)
      {
        comm->docomm[2*idir+1] = comm->buffer_sz[2*idir+1] > 0;
      }
    }
    else if (commflag == PUGH_MINUSFACESCOMM)
    {
      for (idir = 0; idir < GA->extras->dim; idir++)
      {
        comm->docomm[2*idir] = comm->buffer_sz[2*idir] > 0;
      }
    }
    else
    {
      for (idir = 0; idir < GA->extras->dim; idir++)
      {
        if (commflag == PUGH_COMM(idir))
        {
          comm->docomm[2*idir] = comm->buffer_sz[2*idir] > 0;
          comm->docomm[2*idir+1] = comm->buffer_sz[2*idir+1] > 0;
        }
      }
    }

    /* FIXME Add back the check that you have a valid COMM model: Gab */

    /* Handle nsize = 1 type cases. This is only important for one
       processor MPI periodic boundaries */
    for (idir = 0; idir < GA->extras->dim; idir++)
    {
      if (GA->extras->nsize[idir] == 1)
      {
        comm->docomm[2*idir] = 0;
        comm->docomm[2*idir+1] = 0;
      }
    }
  }

#else

  /* get rid of compiler warning about unused parameters */
  pughGH = pughGH;
  comm = comm;
  commflag = commflag;

#endif   /* CCTK_MPI */

  return retval;
}


/*@@
   @routine    PUGH_DisableComm
   @date       Mon Jun 05 2000
   @author     Thomas Radke
   @desc
               This frees the communication buffers of a given comm structure.
   @enddesc
   @history
               Separated from routine PUGH_DisableGArrayDataStorage()
   @endhistory

   @var        GA
   @vdesc      Pointer to grid array structure
   @vtype      pGA *
   @vio        in
   @endvar
   @var        comm
   @vdesc      Pointer to comm structure
   @vtype      pComm *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               1 for communication was disabled
   @endreturndesc
@@*/
static int PUGH_DisableComm(pGA *GA, pComm *comm)
{
#ifdef CCTK_MPI
  int i;                  /* looper */


#ifdef DEBUG_PUGH
  printf (" PUGH_DisableComm: freeing comm buffer for group of %d vars and "
          "first var '%s'\n", comm->n_vars, GA->name);
  fflush (stdout);
#endif

  if (comm->commflag != PUGH_NOCOMM)
  {
    /* free memory for communication buffers: 2 faces per direction */
    for (i = 0; i < 2 * GA->extras->dim; i++)
    {
      if (comm->send_buffer[i])
      {
        free(comm->send_buffer[i]);
        comm->send_buffer[i] = NULL;
      }

      if (comm->recv_buffer[i])
      {
        free(comm->recv_buffer[i]);
        comm->recv_buffer[i] = NULL;
      }

      comm->buffer_sz[i] = 0;
    }

    comm->commflag = PUGH_NOCOMM;
  }
  else
  {
    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                "PUGH_DisableComm: Communication already disabled for "
                "group of %d vars and first var '%s'", comm->n_vars, GA->name);
  }

#else

  /* get rid of compiler warning about unused parameters */
  GA = GA;
  comm = comm;

#endif   /* CCTK_MPI */

  return (1);
}


/*@@
   @routine    PUGH_Sync
   @date       Mon Jun 05 2000
   @author     Thomas Radke
   @desc
               Finally synchronizes a variable or group of variables
               according to a given comm structure.
   @enddesc
   @calls      PUGH_SyncSingleProc

   @var        pughGH
   @vdesc      Pointer to PUGH grid extensions
   @vtype      pGH *
   @vio        in
   @endvar
   @var        comm
   @vdesc      Pointer to comm structure
   @vtype      pComm *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               0 for success
   @endreturndesc
@@*/
static int PUGH_Sync(pGH *pughGH, pComm *comm)
{
#ifdef CCTK_MPI
  int dir;
  pGA *GA;
  MPI_Status mss;
#ifdef PUGH_WITH_DERIVED_DATATYPES
  int i;
  MPI_Request *sr;
#endif
#ifdef COMM_TIMING
  double t1, t2;
#endif
#endif


  /* single-processor case is handled in separate routine */
  if (pughGH->nprocs == 1)
  {
    return (PUGH_SyncSingleProc (pughGH, comm));
  }

#ifdef CCTK_MPI

  /* start the timer for communication time */
  if (pughGH->comm_time >= 0)
  {
    CCTK_TimerStartI (pughGH->comm_time);
  }

  GA = (pGA *) pughGH->variables [comm->first_var][comm->sync_timelevel];

#ifdef PUGH_WITH_DERIVED_DATATYPES
  if (pughGH->commmodel == PUGH_DERIVEDTYPES)
  {
    /* 2 faces, send and receive is the 2 * 2 */
    sr = (MPI_Request *) malloc(comm->n_vars * 2 * 2 * sizeof(MPI_Request));
  }
#endif

#ifdef DEBUG_PUGH
  printf (" PUGH_Sync: syncing group of %d vars with first var '%s'\n",
          comm->n_vars, GA->name);
  fflush (stdout);
#endif

  for (dir = 0; dir < GA->extras->dim; dir ++)
  {

#ifdef COMM_TIMING
    t1 = MPI_Wtime();
#endif

    PostReceiveGA(pughGH, 2*dir, comm);
    PostReceiveGA(pughGH, 2*dir+1, comm);

#ifdef COMM_TIMING
    t2 = MPI_Wtime();
    printf("PR : %f\n",t2-t1);
#endif

    PostSendGA(pughGH, 2*dir, comm);
    PostSendGA(pughGH, 2*dir+1, comm);

#ifdef COMM_TIMING
    t1 = MPI_Wtime();
    printf("PS : %f\n",t1-t2);
#endif

    /* Now comes the big difference between derived types and
       allocated buffers. With derived types, we now have to
       wait on all our recieve AND SEND buffers so we can
       keep on using the send buffers ( as communications are
       in-place). With the allocated we have to wait on each
       recieve, but not on the send, since we don't need the
       send buffer until we pack a send again (above)
    */

    if (pughGH->commmodel == PUGH_ALLOCATEDBUFFERS)
    {
      /* Do a wait any on the receives */
      MPI_Wait(&comm->rreq[2*dir], &mss);
      FinishReceiveGA(pughGH, 2*dir, comm);
      MPI_Wait(&comm->rreq[2*dir+1], &mss);
      FinishReceiveGA(pughGH, 2*dir+1, comm);
    }
#ifdef PUGH_WITH_DERIVED_DATATYPES
    else if (pughGH->commmodel == PUGH_DERIVEDTYPES)
    {
      /* Load up the thing for the waitall */
      for (i = 0; i < comm->n_vars; i++)
      {
        int id = i * 2 * 2;
        pGA *GA = (pGA *) pughGH->variables [i][comm->sync_timelevel];

        if (GA->comm->docomm[2*dir] &&
            GA->storage)
        {
          sr[id] = GA->comm->sreq[2*dir];
          sr[id+1] = GA->comm->rreq[2*dir];
        }
        else
        {
          sr[id] = MPI_REQUEST_NULL;
          sr[id+1] = MPI_REQUEST_NULL;
        }

        if (GA->comm->docomm[2*dir+1] &&
            GA->storage)
        {
          sr[id+2] = GA->comm->sreq[2*dir+1];
          sr[id+3] = GA->comm->rreq[2*dir+1];
        }
        else
        {
          sr[id+2] = MPI_REQUEST_NULL;
          sr[id+3] = MPI_REQUEST_NULL;
        }
      }
      /* Now do a waitall */
      MPI_Waitall(4*comm->n_vars, sr, &mss);
    }
#endif

#ifdef COMM_TIMING
    t2 = MPI_Wtime();
    printf("FR : %f\n",t2-t1);
#endif

  }

#ifdef PUGH_WITH_DERIVED_DATATYPES
  if (pughGH->commmodel == PUGH_DERIVEDTYPES)
  {
    free(sr);
  }
  else
#endif
  {
    /* wait for MPI to finish all outstanding send requests */
    CACTUS_MPI_ERROR (MPI_Waitall (2 * GA->extras->dim, comm->sreq,
                                   comm->sstatus));
  }

  /* get the time spent in communication */
  if (pughGH->comm_time >= 0)
  {
    CCTK_TimerStopI(pughGH->comm_time);
  }

#endif /* CCTK_MPI */

  return (0);
}


/*@@
   @routine    PUGH_SyncSingleProc
   @date       Sun Jan 14 2001
   @author     Thomas Radke
   @desc
               Finally synchronizes a variable or group of variables
               in the single-processor case.
   @enddesc

   @var        pughGH
   @vdesc      Pointer to PUGH grid extensions
   @vtype      pGH *
   @vio        in
   @endvar
   @var        comm
   @vdesc      Pointer to comm structure
   @vtype      pComm *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               0 for success
   @endreturndesc
@@*/
static int PUGH_SyncSingleProc(pGH *pughGH, pComm *comm)
{
  pGA *GA;
  int i, face, dim, copy_bytes, offset_from, offset_to;
  char *data;
  int *istart_from, *iend_from, *iterator_from;
  int *istart_to, *iterator_to;
  int pass, num_passes, pass_total, pass_each, pass_from, pass_to;


  GA = (pGA *) pughGH->variables [comm->first_var][comm->sync_timelevel];

  /* return if no storage assigned */
  if (! GA->storage)
  {
    CCTK_VWarn(2, __LINE__, __FILE__, CCTK_THORNSTRING,
               "Trying to synchronize variable '%s' with no storage", GA->name);
    return (0);
  }

  /* since we need two iterators here we need to allocate one */
  iterator_from = GA->extras->iterator;
  iterator_to = (int *) malloc (GA->extras->dim * sizeof (int));

  /* loop over all faces */
  for (face = 0; face < 2 * GA->extras->dim; face++)
  {

    /* there's only something to do for periodic boundary conditions */
    if (GA->connectivity->perme[face / 2])
    {
      /* if the domain is smaller than the number of ghost zones, then
         we have to copy in several passes */
      /* FIXME: these passes do not each need to copy the whole
         regaion, they could copy a part of it only */
      /* total number of points to copy */
      pass_total = GA->extras->nghostzones[face>>1];
      /* points copied in one pass */
      pass_each = (GA->extras->lnsize[face>>1]
                   - 2 * GA->extras->nghostzones[face>>1]);
      /* necessary number of passes */
      num_passes = (pass_total + pass_each - 1) / pass_each;
      for (pass = 0; pass < num_passes; pass++)
      {
        /* calculate what needs to be done in this pass */
        /* FIXME: these values are currently ignored */
        pass_from = pass * pass_each;
        pass_to = pass_from + pass_each - 1;
        if (pass_to > pass_total) pass_to = pass_total;

        /* get the index ranges for the nested loops */
        istart_from   = GA->extras->overlap[GA->stagger[face>>1]][0][face];
        iend_from     = GA->extras->overlap[GA->stagger[face>>1]][1][face];
        if (face & 1)
        {
          istart_to   = GA->extras->ghosts[GA->stagger[face>>1]][0][face - 1];
        }
        else
        {
          istart_to   = GA->extras->ghosts[GA->stagger[face>>1]][0][face + 1];
        }

        /* set iterators to the start vector */
        for(dim = 0; dim < GA->extras->dim; dim++)
        {
          iterator_from[dim] = istart_from[dim];
          iterator_to[dim] = istart_to[dim];
        }

        /* get the number of bytes to copy in the lowest dimension */
        copy_bytes = (iend_from[0] - istart_from[0]) * GA->varsize;

        /* now do the nested loops starting with the innermost */
        dim = 1;
        while (1)
        {
          /* check for end of current loop */
          if (GA->extras->dim > 1 && iterator_from[dim] >= iend_from[dim])
          {
            /* increment outermost loopers */
            for (dim++; dim < GA->extras->dim; dim++)
            {
              iterator_from[dim]++;
              iterator_to[dim]++;
              if (iterator_from[dim] < iend_from[dim])
              {
                break;
              }
            }

            /* done if beyond outermost loop */
            if (dim >= GA->extras->dim)
            {
              break;
            }

            /* reset innermost loopers */
            for (dim--; dim > 0; dim--)
            {
              iterator_from[dim] = istart_from[dim];
              iterator_to[dim] = istart_to[dim];
            }
            dim = 1;
          }

          /* get the linear offsets */
          offset_from = istart_from[0];
          offset_to   = istart_to[0];
          for (i = 1; i < GA->extras->dim; i++)
          {
            offset_from += iterator_from[i] * GA->extras->hyper_volume[i];
            offset_to   += iterator_to[i]   * GA->extras->hyper_volume[i];
          }
          offset_from *= GA->varsize;
          offset_to   *= GA->varsize;

          /* copy the data for all the variables in the comm structure */
          for (i = comm->first_var; i < comm->first_var + comm->n_vars; i++)
          {
            data = ((pGA *) pughGH->variables[i][comm->sync_timelevel])->data;
            memcpy (data + offset_to, data + offset_from, copy_bytes);
          }

          /* increment current loopers */
          if (GA->extras->dim > 1)
          {
            /* increment current looper */
            iterator_from[dim]++;
            iterator_to[dim]++;
          }
          else
          {
            /* exit loop if array dim is 1D */
            break;
          }

        } /* end of nested loops over all dimensions */
      } /* for pass */
    } /* if periodic */
  } /* for face */

  /* free the allocated iterator */
  free (iterator_to);

  return (0);
}