aboutsummaryrefslogtreecommitdiff
path: root/src/Comm.c
blob: d4614b5c1ca50f8f76b1f0537d5ba54f348aac77 (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
 /*@@
   @file      Comm.c
   @date      Thu Feb  4 11:34:29 1999
   @author    Tom Goodale
   @desc 
   Pugh communication functions
   @enddesc 
   @version $Header$
 @@*/

/*#define DEBUG_PUGH 1*/

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

#include "cctk.h"

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

static 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(pGH *pughGH,
                            pComm *comm);
static int PUGH_SyncGArrayGroup(pGH *pughGH,
                                int first_var);
static int PUGH_Sync(pGH *pughGH,
                     pComm *comm);


/*@@
  @routine PUGH_SyncGroup
  @author  Thomas Radke
  @date    30 Mar 1999
  @desc
  Synchronizes all variables in the group indicated by groupname.
  Only groups of type GROUP_ARRAY and GROUP_GF can be synchronized.
  @enddesc
  @calls   CCTK_DecomposeName CCTK_GroupIndex CCTK_GroupData CCTK_WARN
           PUGH_SyncGArrayGroup
  @history
 
  @endhistory 
  @var     GH
  @vdesc   Pointer to CCTK grid hierarchy
  @vtype   cGH
  @vio     in
  @endvar 
  @var     groupname
  @vdesc   name of the group to be synchronized
  @vtype   const char *
  @vio     in
  @endvar 
  @@*/

int PUGH_SyncGroup(cGH *GH, const char *groupname)
{
  cGroup pgroup;           /* group information */
  int group;               /* group index */
  int rc;                  /* return code */

#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);
  CCTK_GroupData(group, &pgroup);

  if (pgroup.grouptype == CCTK_SCALAR)
  {
    rc = 0;
    CCTK_WARN(4, "Synchronising a scalar in PUGH");
  }
  else if (pgroup.grouptype == CCTK_GF || pgroup.grouptype == CCTK_ARRAY)
  {
    rc = PUGH_SyncGArrayGroup(PUGH_pGH(GH), CCTK_FirstVarIndexI(group));
  }
  else
  {
    CCTK_WARN(1, "Unknown group type in PUGH_SyncGroup");
    rc = 0;
  }

  return (rc);
}


/*@@
  @routine PUGH_EnableGroupComm
  @author  Thomas Radke
  @date    30 Mar 1999
  @desc
  Enables communication for all variables in the group indicated by groupname.
  @enddesc
  @calls   CCTK_DecomposeName CCTK_GroupIndex CCTK_GroupData CCTK_WARN
  @history
 
  @endhistory 
  @var     GH
  @vdesc   Pointer to CCTK grid hierarchy
  @vtype   cGH
  @vio     in
  @endvar 
  @var     groupname
  @vdesc   name of the group to be synchronized
  @vtype   const char *
  @vio     in
  @endvar 
  @@*/

int PUGH_EnableGroupComm(cGH *GH, const char *groupname)
{
  int group;               /* group index */
  cGroup pgroup;           /* group information */
  int rc;                  /* return code */

#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);

  if (pgroup.grouptype == CCTK_SCALAR)
  {
    rc = 1;
  }
  else if (pgroup.grouptype == CCTK_GF || pgroup.grouptype == CCTK_ARRAY)
  {
    rc = PUGH_EnableGArrayGroupComm(PUGH_pGH(GH),
                                    CCTK_FirstVarIndexI(group),
                                    PUGH_ALLCOMM);
  }
  else
  {
    CCTK_WARN(1, "Unknown group type in PUGH_EnableGroupComm");
    rc = 0;
  }

  return (rc);
}


/*@@
  @routine PUGH_DisableGroupComm
  @author  Thomas Radke
  @date    30 Mar 1999
  @desc
  Disables communication for all variables in the group indicated by groupname.
  @enddesc
  @calls   CCTK_DecomposeName CCTK_GroupIndex CCTK_GroupData CCTK_WARN
  @history
 
  @endhistory 
  @var     GH
  @vdesc   Pointer to CCTK grid hierarchy
  @vtype   cGH
  @vio     in
  @endvar 
  @var     groupname
  @vdesc   name of the group to be synchronized
  @vtype   const char *
  @vio     in
  @endvar 
  @@*/

int PUGH_DisableGroupComm(cGH *GH, const char *groupname)
{
  int group;               /* group index */
  cGroup pgroup;           /* pointer to group information */
  int rc;                  /* return code */

  pGH *pughGH;
  int var;

#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);

  if (pgroup.grouptype == CCTK_SCALAR)
  {
    rc = 1;
  }
  else if (pgroup.grouptype == CCTK_GF || pgroup.grouptype == CCTK_ARRAY)
  {
    pughGH=PUGH_pGH(GH);
    var = CCTK_FirstVarIndexI(group);

    /* FIXME:  workaround.  This one is really bad ! */
    rc = PUGH_DisableGArrayGroupComm(pughGH, var,(((pGA ***)pughGH->variables)[var][0])->groupcomm);
  }
  else
  {
    CCTK_WARN(1, "Unknown group type in PUGH_DisableGroupComm");
    rc = 0;
  }

  return (rc);
}


 /*@@
   @routine    PUGH_EnableGArrayComm
   @date       Mon Jun 05 2000
   @author     Thomas Radke
   @desc
               Enables communication for a single array.
   @enddesc
   @history
   @endhistory
@@*/
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
   @history
   @endhistory
@@*/
int PUGH_DisableGArrayComm(pGA *GA)
{

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

  return (PUGH_DisableComm((pGH *) GA->parent, 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
   @history
   @endhistory
@@*/
int PUGH_SyncGArrayGroup(pGH *pughGH,
                         int first_var)
{
  pGA *firstGA;


  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

  return (PUGH_Sync(pughGH, firstGA->groupcomm));
}


 /*@@
   @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
   @history
   @endhistory
@@*/
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));
}


int PUGH_Barrier(cGH *GH)
{
#ifdef CCTK_MPI
  CACTUS_MPI_ERROR(MPI_Barrier(PUGH_pGH(GH)->PUGH_COMM_WORLD));
#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
   @history
   @endhistory
@@*/
static int PUGH_EnableGArrayGroupComm(pGH *pughGH,
                                      int first_var,
                                      int commflag)
{
  pGA *GA;                 /* first variable in group */


  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_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
   @history
   @endhistory
@@*/
int PUGH_DisableGArrayGroupComm(pGH *pughGH,
                                int first_var,
                                pComm *groupcomm)
{
#ifdef DEBUG_PUGH
  pGA *GA;                 /* first variable in group */

  GA = (pGA *) pughGH->variables[first_var][0];
  printf(" PUGH_DisableGArrayGroupComm: request for group "
         "with first var '%s'\n", GA->name);
  fflush(stdout);
#endif

  return (PUGH_DisableComm(pughGH, groupcomm));
}


 /*@@
   @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
@@*/
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++)
    {
      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--)
            {
              if (comm->send_buffer[i])
              {
                free(comm->send_buffer[i]);
              }
              if (comm->recv_buffer[i])
              {
                free(comm->recv_buffer[i]);
              }
              comm->buffer_sz[i] = 0;
            }

            CCTK_WARN(1, "Out of memory !");
            retval = -1;
            break;
          }
        }
        else
        {
          comm->buffer_sz[i]   = 0;
          comm->send_buffer[i] = NULL;
          comm->recv_buffer[i] = NULL;
        }
      }
      else               /* no neighbor -> no comm buffers */
      {
        comm->buffer_sz[i]   = 0;
        comm->send_buffer[i] = NULL;
        comm->recv_buffer[i] = NULL;
      }
    }
  }

  /* 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;
      }
    }
  }
#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
@@*/
static int PUGH_DisableComm(pGH *pughGH,
                            pComm *comm)
{
#ifdef CCTK_MPI
  int i;                  /* looper */
  pGA *GA;                /* GA structure the comm structure belongs to */


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

#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,
                "Communication already disabled for group of %d vars "
                "and first var '%s'", comm->n_vars, GA->name);
  }
#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
   @history
   @endhistory
@@*/
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


#if 0
  /* start the timer for communication time */
  CCTK_TimerStartI(pughGH->comm_time);
#endif

#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

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

#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
  free(sr);
#endif

#if 0
  /* get the time spent in communication */
  CCTK_TimerStopI(pughGH->comm_time);
#endif

#endif            /* CCTK_MPI */

  return (0);
}