Пример #1
0
int
main (int argc, char** argv)
{
  QMP_status_t status;
  int this_node;
  QMP_thread_level_t req, prv;

  /* Start QMP */
  req = QMP_THREAD_SINGLE;
  status = QMP_init_msg_passing (&argc, &argv, req, &prv);

  if (status != QMP_SUCCESS) {
    QMP_error ("QMP_init failed: %s\n", QMP_error_string(status));
    QMP_abort(1);
  }

  /* Get my logical node number */
  this_node = QMP_get_node_number();

  /* Print the result */
  printf("%04d",this_node);

  /* Quit */
  QMP_finalize_msg_passing ();

  return 0;
}
Пример #2
0
/*--------------------------------------------------------------------*/
static void setup_qmp_grid(){
  int ndim = 4;
  int len[4];
  int ndim2, i;
  const int *nsquares2;

  len[0] = nx; len[1] = ny; len[2] = nz; len[3] = nt;

  if(mynode()==0){
    printf("qmp_grid,");
    printf("\n");
  }

  ndim2 = QMP_get_allocated_number_of_dimensions();
  nsquares2 = QMP_get_allocated_dimensions();

  /* If the dimensions are not already allocated, use the
     node_geometry request.  Otherwise a hardware or command line
     specification trumps the parameter input. */
#ifdef FIX_NODE_GEOM
  if(ndim2 == 0){
    ndim2 = 4;
    nsquares2 = node_geometry;
  }
  else{
    node0_printf("setup_qmp_grid: Preallocated machine geometry overrides request\n");
  }
#endif

  if(mynode()==0){
    printf("Using machine geometry: ");
    for(i=0; i<ndim; i++){
      printf("%d ",nsquares2[i]);
      if(i < ndim-1)printf("X ");
    }
    printf("\n");
  }

  /* In principle, we could now rotate coordinate axes */
  /* Save this for a future upgrade */

  set_qmp_layout_grid(nsquares2, ndim2);

  ndim2 = QMP_get_logical_number_of_dimensions();
  nsquares2 = QMP_get_logical_dimensions();

  for(i=0; i<ndim; i++) {
    if(i<ndim2) nsquares[i] = nsquares2[i];
    else nsquares[i] = 1;
  }

  for(i=0; i<ndim; i++) {
    if(len[i]%nsquares[i] != 0) {
      node0_printf("LATTICE SIZE DOESN'T FIT GRID\n");
      QMP_abort(0);
    }
    squaresize[i] = len[i]/nsquares[i];
  }
}
Пример #3
0
void comm_abort(int status)
{
  QMP_abort(status);
}
Пример #4
0
void comm_exit(int ret)
{
  if (ret) QMP_abort(ret);
}
Пример #5
0
int
main (int argc, char** argv)
{
  int             i, nc;
  QMP_status_t      status;
  int       **smem, **rmem;
  QMP_msgmem_t    *recvmem;
  QMP_msghandle_t *recvh;
  QMP_msgmem_t    *sendmem;
  QMP_msghandle_t *sendh;
  struct perf_argv pargv;
  QMP_thread_level_t req, prv;

  /** 
   * Simple point to point topology 
   */
  int dims[4] = {2,2,2,2};
  int ndims = 1;

  //if(QMP_get_node_number()==0)
  //printf("starting init\n"); fflush(stdout);
  req = QMP_THREAD_SINGLE;
  status = QMP_init_msg_passing (&argc, &argv, req, &prv);
  if (status != QMP_SUCCESS) {
    fprintf (stderr, "QMP_init failed\n");
    return -1;
  }
  if(QMP_get_node_number()==0)
    printf("finished init\n"); fflush(stdout);

  if (parse_options (argc, argv, &pargv) == -1) {
    if(QMP_get_node_number()==0)
      usage (argv[0]);
    exit (1);
  }

  {
    int maxdims = 4;
    int k=0;
    int nodes = QMP_get_number_of_nodes();
    ndims = 0;
    while( (nodes&1) == 0 ) {
      if(ndims<maxdims) ndims++;
      else {
	dims[k] *= 2;
	k++;
	if(k>=maxdims) k = 0;
      }
      nodes /= 2;
    }
    if(nodes != 1) {
      QMP_error("invalid number of nodes %i", QMP_get_number_of_nodes());
      QMP_error(" must power of 2");
      QMP_abort(1);
    }
    pargv.ndims = ndims;
  }

  status = QMP_declare_logical_topology (dims, ndims);
  if (status != QMP_SUCCESS) {
    fprintf (stderr, "Cannot declare logical grid\n");
    return -1;
  }

  /* do a broadcast of parameter */
  if (QMP_broadcast (&pargv, sizeof (pargv)) != QMP_SUCCESS) {
    QMP_printf ("Broadcast parameter failed\n");
    exit (1);
  }

  {
    int k=1;
    const int *lc = QMP_get_logical_coordinates();
    for(i=0; i<ndims; i++) k += lc[i];
    pargv.sender = k&1;
  }

  QMP_printf("%s options: num_channels[%d] verify[%d] option[%d] datasize[%d] numloops[%d] sender[%d] strided_send[%i] strided_recv[%i] strided_array_send[%i] ",
	     argv[0], pargv.num_channels, pargv.verify, 
	     pargv.option, pargv.size, pargv.loops, pargv.sender,
	     strided_send, strided_recv, strided_array_send);
  fflush(stdout);


  /**
   * Create memory
   */
  nc = pargv.num_channels;
  smem = (int **)malloc(nc*sizeof (int *));
  rmem = (int **)malloc(nc*sizeof (int *));
  sendmem = (QMP_msgmem_t *)malloc(ndims*nc*sizeof (QMP_msgmem_t));
  recvmem = (QMP_msgmem_t *)malloc(ndims*nc*sizeof (QMP_msgmem_t));
  sendh = (QMP_msghandle_t *)malloc(nc*sizeof (QMP_msghandle_t));
  recvh = (QMP_msghandle_t *)malloc(nc*sizeof (QMP_msghandle_t));

  QMP_barrier();
  if(QMP_get_node_number()==0) printf("\n"); fflush(stdout);
  if(pargv.option & TEST_SIMUL) {
    int opts = pargv.option;
    pargv.option = TEST_SIMUL;
    if(QMP_get_node_number()==0)
      QMP_printf("starting simultaneous sends"); fflush(stdout);
    for(i=pargv.minsize; i<=pargv.maxsize; i*=pargv.facsize) {
      pargv.size = i;
      create_msgs(smem, rmem, sendmem, recvmem, sendh, recvh, ndims, nc, i, &pargv);
      test_simultaneous_send (smem, rmem, sendh, recvh, &pargv);
      check_mem(rmem, ndims, nc, i);
      free_msgs(smem, rmem, sendmem, recvmem, sendh, recvh, ndims, nc);
    }
    if(QMP_get_node_number()==0)
      QMP_printf("finished simultaneous sends\n"); fflush(stdout);
    pargv.option = opts;
  }

  if(pargv.option & TEST_PINGPONG) {
    int opts = pargv.option;
    pargv.option = TEST_PINGPONG;
    if(QMP_get_node_number()==0)
      QMP_printf("starting ping pong sends"); fflush(stdout);
    for(i=pargv.minsize; i<=pargv.maxsize; i*=pargv.facsize) {
      pargv.size = i;
      create_msgs(smem, rmem, sendmem, recvmem, sendh, recvh, ndims, nc, i, &pargv);
      if(pargv.verify)
	test_pingpong_verify(smem, rmem, sendh, recvh, &pargv);
      else
	test_pingpong(smem, rmem, sendh, recvh, &pargv);
      check_mem(rmem, ndims, nc, i);
      free_msgs(smem, rmem, sendmem, recvmem, sendh, recvh, ndims, nc);
    }
    if(QMP_get_node_number()==0)
      QMP_printf("finished ping pong sends\n"); fflush(stdout);
    pargv.option = opts;
  }

  if(pargv.option & TEST_ONEWAY) {
    int opts = pargv.option;
    pargv.option = TEST_ONEWAY;
    if(QMP_get_node_number()==0)
      QMP_printf("starting one way sends"); fflush(stdout);
    for(i=pargv.minsize; i<=pargv.maxsize; i*=pargv.facsize) {
      pargv.size = i;
      create_msgs(smem, rmem, sendmem, recvmem, sendh, recvh, ndims, nc, i, &pargv);
      test_oneway (smem, rmem, sendh, recvh, &pargv);
      if(!pargv.sender) check_mem(rmem, ndims, nc, i);
      free_msgs(smem, rmem, sendmem, recvmem, sendh, recvh, ndims, nc);
    }
    if(QMP_get_node_number()==0)
      QMP_printf("finished one way sends"); fflush(stdout);
    pargv.option = opts;
  }


  /**
   * Free memory 
   */
  free (smem);
  free (rmem);

  free (sendh);
  free (recvh);

  free (sendmem);
  free (recvmem);

  QMP_finalize_msg_passing ();

  return 0;
}
Пример #6
0
void
create_msgs(int **smem, int **rmem,
	    QMP_msgmem_t *sendmem, QMP_msgmem_t *recvmem,
	    QMP_msghandle_t *sendh, QMP_msghandle_t *recvh,
	    int ndims, int nc, int size, struct perf_argv *pargv)
{
  int i, j, n;

  for (i = 0; i < nc; i++) {

    if(strided_array_send) {
      void *base[NAMAX];
      size_t bsize[NAMAX];
      int nblocks[NAMAX];
      ptrdiff_t stride[NAMAX];
      int tsize, skip;
      int na, k, bs, nb, nbt, ab, as, st;

      bs = strided_array_send;
      nbt = size/bs;
      na = sqrt(nbt);
      if(na<2) na = nbt;
      if(na>NAMAX) na = NAMAX;
      nb = nbt/na;
      st = 2*bs;
      skip = 3*bs;
      ab = bs*nb;
      as = st*nb+skip;

      tsize = 0;
      for(k=0; k<na; k++) {
	bsize[k] = bs*sizeof(int);
	stride[k] = st*sizeof(int);
	nblocks[k] = nb;
	if(k==na-1) nblocks[k] = nbt - nb*(na-1);
	tsize += skip + st * nblocks[k];
      }

      smem[i] = (int *)malloc(ndims*tsize*sizeof(int));
      for(n=0; n<ndims; n++) {
	for (j = 0; j < tsize; j++) { smem[i][n*tsize+j] = 0; }
	for (j = 0; j < size; j++) {
	  int ai, ak;
	  ak = j/ab;
	  if(ak>=na) ak = na-1;
	  ai = j-(ab*ak);
	  k = (as*ak) + st*(ai/bs) + (ai%bs);
	  smem[i][n*tsize+k] = i+j+1;
	}
	for(k=0; k<na; k++) {
	  base[k] = (void *)&(smem[i][n*tsize+as*k]);
	}
	sendmem[n*nc+i] =
	  QMP_declare_strided_array_msgmem(base, bsize, nblocks, stride, na);
	if(!sendmem[n*nc+i]) {
	  QMP_printf("error in declare strided msgmem\n");
	  QMP_abort(1);
	}
      }
    } else
    if(strided_send) {
      int tsize, bsize, stride, nblocks;

      bsize = strided_send;
      stride = 2*bsize;
      nblocks = size/bsize;
      tsize = stride * nblocks;

      smem[i] = (int *)malloc(ndims*tsize*sizeof (int));
      for(n=0; n<ndims; n++) {
	for (j = 0; j < tsize; j++) { smem[i][n*tsize+j] = 0; }
	for (j = 0; j < size; j++) {
	  int k = stride*(j/bsize) + (j%bsize);
	  smem[i][n*tsize+k] = i+j+1;
	}
	sendmem[n*nc+i] = QMP_declare_strided_msgmem(smem[i]+(n*tsize),
						     bsize*sizeof(int),
						     nblocks,
						     stride*sizeof(int));
	if(!sendmem[n*nc+i]) {
	  QMP_printf("error in declare strided msgmem\n");
	  QMP_abort(1);
	}
      }

    } else {

      smem[i] = (int *)malloc(ndims*size*sizeof(int));
      for(n=0; n<ndims; n++) {
	for (j = 0; j < size; j++) {
	  smem[i][n*size+j] = i+j+1;
	}
	sendmem[n*nc+i] = QMP_declare_msgmem(smem[i]+(n*size),
					     size*sizeof(int));
	if(!sendmem[n*nc+i]) {
	  QMP_printf("error in declare msgmem\n");
	  QMP_abort(1);
	}
      }

    }

    if(strided_recv) {
      int tsize, bsize, stride, nblocks;

      bsize = strided_recv;
      stride = 2*bsize;
      nblocks = size/bsize;
      tsize = stride * nblocks;

      rmem[i] = (int *)malloc(ndims*tsize*sizeof (int));
      for(n=0; n<ndims; n++) {
	for (j = 0; j < tsize; j++) {
	  rmem[i][n*tsize+j] = 0;
	}
	recvmem[n*nc+i] = QMP_declare_strided_msgmem(rmem[i]+(n*tsize),
						     bsize*sizeof(int),
						     nblocks,
						     stride*sizeof(int));
	if(!recvmem[n*nc+i]) {
	  QMP_printf("error in declare strided msgmem\n");
	  QMP_abort(1);
	}
      }

    } else {

      rmem[i] = (int *)malloc(ndims*size*sizeof (int));
      for(n=0; n<ndims; n++) {
	for (j = 0; j < size; j++) {
	  rmem[i][n*size+j] = 0;
	}
	recvmem[n*nc+i] = QMP_declare_msgmem (rmem[i]+(n*size),
					      size*sizeof(int));
	if(!recvmem[n*nc+i]) {
	  QMP_printf("error in declare msgmem\n");
	  QMP_abort(1);
	}
      }

    }

    if(ndims>0) {  // always use
      QMP_msghandle_t *tsend, *trecv;
      tsend = (QMP_msghandle_t *)malloc(ndims*sizeof(QMP_msghandle_t));
      trecv = (QMP_msghandle_t *)malloc(ndims*sizeof(QMP_msghandle_t));
      for(n=0; n<ndims; n++) {
	trecv[n] = QMP_declare_receive_relative (recvmem[n*nc+i], n, 1, 0);
	if (!trecv[n]) {
	  QMP_printf ("Recv Handle Error: %s\n", QMP_get_error_string(0));
	  exit (1);
	}
	tsend[n] = QMP_declare_send_relative (sendmem[n*nc+i], n, -1, 0);
	if (!tsend[n]) {
	  QMP_printf ("Send Handle Error: %s\n", QMP_get_error_string(0));
	  exit (1);
	}
      }
      if(pargv->option & TEST_PINGPONG) {
	if(pargv->sender) {
	  sendh[i] = QMP_declare_send_recv_pairs(tsend, ndims);
	  recvh[i] = QMP_declare_send_recv_pairs(trecv, ndims);
	} else {
	  recvh[i] = QMP_declare_send_recv_pairs(trecv, ndims);
	  sendh[i] = QMP_declare_send_recv_pairs(tsend, ndims);
	}
      } else {
	recvh[i] = QMP_declare_multiple(trecv, ndims);
	sendh[i] = QMP_declare_multiple(tsend, ndims);
      }
      if (!recvh[i]) {
	QMP_printf ("Recv Handle Error: %s\n", QMP_get_error_string(0));
	exit (1);
      }
      if (!sendh[i]) {
	QMP_printf ("Send Handle Error: %s\n", QMP_get_error_string(0));
	exit (1);
      }
      free(tsend);
      free(trecv);
    } else {
      recvh[i] = QMP_declare_receive_relative (recvmem[i], 0, 1, 0);
      if (!recvh[i]) {
	QMP_printf ("Recv Handle Error: %s\n", QMP_get_error_string(0));
	exit (1);
      }
      sendh[i] = QMP_declare_send_relative (sendmem[i], 0, -1, 0);
      if (!sendh[i]) {
	QMP_printf ("Send Handle Error: %s\n", QMP_get_error_string(0));
	exit (1);
      }
    }

  }
}
void make_shift_tables(int bound[2][4][4], halfspinor_array* chi1,
                       halfspinor_array* chi2,

                       halfspinor_array* recv_bufs[2][4],
                       halfspinor_array* send_bufs[2][4],

                       void (*QDP_getSiteCoords)(int coord[], int node, int linearsite),
                       int (*QDP_getLinearSiteIndex)(const int coord[]),

                       int (*QDP_getNodeNumber)(const int coord[]))
{
    volatile int dir,i;
    const int my_node = QMP_get_node_number();
    int coord[4];
    int gcoord[4];
    int gcoord2[4];

    int linear;
    int **shift_table;
    int x,y,z,t;
    int *subgrid_size = getSubgridSize();
    int mu;
    int offset;

    int cb;
    const int *node_coord  = QMP_get_logical_coordinates();
    int p;
    int site, index;

    InvTab4 *xinvtab;
    InvTab4 *invtab;

    int qdp_index;
    int my_index;
    int num;
    int offsite_found;

    /* Setup the subgrid volume for ever after */
    subgrid_vol = 1;
    for(i=0; i < getNumDim(); ++i) {
        subgrid_vol *= getSubgridSize()[i];
    }

    /* Get the checkerboard size for ever after */
    subgrid_vol_cb = subgrid_vol / 2;

    /* Now I want to build the site table */
    /* I want it cache line aligned? */
    xsite_table = (int *)malloc(sizeof(int)*subgrid_vol+63L);
    if(xsite_table == 0x0 ) {
        QMP_error("Couldnt allocate site table");
        QMP_abort(1);
    }

    site_table = (int *)((((ptrdiff_t)(xsite_table))+63L)&(-64L));

    xinvtab = (InvTab4 *)malloc(sizeof(InvTab4)*subgrid_vol + 63L);
    if(xinvtab == 0x0 ) {
        QMP_error("Couldnt allocate site table");
        QMP_abort(1);
    }
    invtab = (InvTab4 *)((((ptrdiff_t)(xinvtab))+63L)&(-64L));

    /* Inversity of functions check:
       Check that myLinearSiteIndex3D is in fact the inverse
       of mySiteCoords3D, and that QDP_getSiteCoords is the
       inverse of QDP_linearSiteIndex()
    */
    for(p=0; p < 2; p++) {
        for(site=0; site < subgrid_vol_cb; site++) {

            /* Linear site index */
            my_index = site + subgrid_vol_cb*p;
            QDP_getSiteCoords(gcoord, my_node, my_index);
            linear=QDP_getLinearSiteIndex(gcoord);

            if( linear != my_index ) {
                printf("P%d cb=%d site=%d : QDP_getSiteCoords not inverse of QDP_getLinearSiteIndex(): my_index=%d linear=%d\n", my_node, p,site, my_index,linear);
            }

            mySiteCoords4D(gcoord, my_node, my_index);
            linear=myLinearSiteIndex4D(gcoord);

            if( linear != my_index ) {
                printf("P%d cb=%d site=%d : mySiteCoords3D not inverse of myLinearSiteIndex3D(): my_index=%d linear=%d\n", my_node, p,site, my_index,linear);
            }
        }
    }


    /* Loop through sites - you can choose your path below */
    /* This is a checkerboarded order which is identical hopefully
       to QDP++'s rb2 subset when QDP++ is in a CB2 layout */
    for(p=0; p < 2; p++) {
        for(t=0; t < subgrid_size[3]; t++) {
            for(z=0; z < subgrid_size[2]; z++) {
                for(y=0; y < subgrid_size[1]; y++) {
                    for(x=0; x < subgrid_size[0]/2; x++) {

                        coord[0] = 2*x + p;
                        coord[1] = y;
                        coord[2] = z;
                        coord[3] = t;

                        /* Make global */
                        for(i=0; i < 4; i++) {
                            coord[i] += subgrid_size[i]*node_coord[i];
                        }

                        /* Index of coordinate -- NB this is not lexicographic
                           but takes into account checkerboarding in QDP++ */
                        qdp_index = QDP_getLinearSiteIndex(coord);

                        /* Index of coordinate in my layout. -- NB this is not lexicographic
                           but takes into account my 3D checkerbaording */
                        my_index = myLinearSiteIndex4D(coord);
                        site_table[my_index] = qdp_index;

                        cb=parity(coord);
                        linear = my_index%subgrid_vol_cb;

                        invtab[qdp_index].cb=cb;
                        invtab[qdp_index].linearcb=linear;
                    }
                }
            }
        }
    }

    /* Site table transitivity check:
       for each site, convert to index in cb3d, convert to qdp index
       convert qdp_index to coordinate
       convert coordinate to back index in cb3d
       Check that your cb3d at the end is the same as you
       started with */
    for(p=0; p < 2; p++) {
        for(site=0; site < subgrid_vol_cb; site++) {

            /* My local index */
            my_index = site + subgrid_vol_cb*p;

            /* Convert to QDP index */
            qdp_index = site_table[ my_index ];

            /* Switch QDP index to coordinates */
            QDP_getSiteCoords(gcoord, my_node,qdp_index);

            /* Convert back to cb3d index */
            linear = myLinearSiteIndex4D(gcoord);

            /* Check new cb,cbsite index matches the old cb index */
            if (linear != my_index) {
                printf("P%d The Circle is broken. My index=%d qdp_index=%d coords=%d,%d,%d,%d linear(=my_index?)=%d\n", my_node, my_index, qdp_index, gcoord[0],gcoord[1],gcoord[2],gcoord[3],linear);
            }
        }
    }


    /* Consistency check 2: Test mySiteCoords 3D
       for all 3d cb,cb3index convert to
       cb3d linear index (my_index)
       convert to qdp_index (lookup in site table)

       Now convert qdp_index and my_index both to
       coordinates. They should produce the same coordinates
    */
    for(p=0; p < 2; p++) {
        for(site=0; site < subgrid_vol_cb; site++) {

            /* My local index */
            my_index = site + subgrid_vol_cb*p;
            mySiteCoords4D(gcoord, my_node, my_index);

            qdp_index = site_table[ my_index ];
            QDP_getSiteCoords(gcoord2, my_node,qdp_index);

            for(mu=0 ; mu < 4; mu++) {
                if( gcoord2[mu] != gcoord[mu] ) {
                    printf("P%d: my_index=%d qdp_index=%d mySiteCoords=(%d,%d,%d,%d) QDPsiteCoords=(%d,%d,%d,%d)\n", my_node, my_index, qdp_index, gcoord[0], gcoord[1], gcoord[2], gcoord[3], gcoord2[0], gcoord2[1], gcoord2[2], gcoord2[3]);
                    continue;
                }
            }

        }
    }

    /* Allocate the shift table */
    /* The structure is as follows: There are 4 shift tables in order:

      [ Table 1 | Table 2 | Table 3 | Table 4 ]
      Table 1: decomp_scatter_index[mu][site]
      Table 2: decomp_hvv_scatter_index[mu][site]
      Table 3: recons_mvv_gather_index[mu][site]
      Table 4: recons_gather_index[mu][site]

    */

    /* This 4 is for the 4 tables: Table 1-4*/
    if ((shift_table = (int **)malloc(4*sizeof(int*))) == 0 ) {
        QMP_error("init_wnxtsu3dslash: could not initialize shift_table");
        QMP_abort(1);

    }

    for(i=0; i < 4; i++) {
        /* This 4 is for the 4 comms dierctions: */
        if ((shift_table[i] = (int *)malloc(4*subgrid_vol*sizeof(int))) == 0) {
            QMP_error("init_wnxtsu3dslash: could not initialize shift_table");
            QMP_abort(1);
        }
    }


    /* Initialize the boundary counters */
    for(cb=0; cb < 2; cb++) {
        for(dir=0; dir < 4; dir++) {
            bound[cb][0][dir] = 0;
            bound[cb][1][dir] = 0;
            bound[cb][2][dir] = 0;
            bound[cb][3][dir] = 0;
        }
    }


    for(cb=0; cb < 2; cb++) {
        for(site=0; site < subgrid_vol_cb; ++site) {

            index = cb*subgrid_vol_cb + site;

            /* Fetch site from site table */
            qdp_index = site_table[index];

            /* Get its coords */
            QDP_getSiteCoords(coord, my_node, qdp_index);

            /* Loop over directions building up shift tables */
            for(dir=0; dir < 4; dir++) {

                int fcoord[4], bcoord[4];
                int fnode, bnode;
                int blinear, flinear;

                /* Backwards displacement*/
                offs(bcoord, coord, dir, -1);
                bnode   = QDP_getNodeNumber(bcoord);
                blinear = QDP_getLinearSiteIndex(bcoord);

                /* Forward displacement */
                offs(fcoord, coord, dir, +1);
                fnode   = QDP_getNodeNumber(fcoord);
                flinear = QDP_getLinearSiteIndex(fcoord);

                /* Scatter:  decomp_{plus,minus} */
                /* Operation: a^F(shift(x,type=0),dir) <- decomp(psi(x),dir) */
                /* Send backwards - also called a receive from forward */
                if (bnode != my_node) {
                    /* Offnode */
                    /* Append to Tail 1, increase boundary count */
                    /* This is the correct code */
                    shift_table[DECOMP_SCATTER][dir+4*index]
                        = subgrid_vol_cb + bound[1-cb][DECOMP_SCATTER][dir];

                    bound[1-cb][DECOMP_SCATTER][dir]++;

                }
                else {
                    /* On node. Note the linear part of its (cb3, linear) bit,
                       using a reverse lookup */
                    shift_table[DECOMP_SCATTER][dir+4*index] =
                        invtab[blinear].linearcb;
                }


                /* Scatter:  decomp_hvv_{plus,minus} */
                /* Operation:  a^B(shift(x,type=1),dir) <- U^dag(x,dir)*decomp(psi(x),dir) */
                /* Send forwards - also called a receive from backward */
                if (fnode != my_node) {
                    /* Offnode */
                    /* Append to Tail 1, increase boundary count */
                    shift_table[DECOMP_HVV_SCATTER][dir+4*index]
                        = subgrid_vol_cb + bound[1-cb][DECOMP_HVV_SCATTER][dir];

                    bound[1-cb][DECOMP_HVV_SCATTER][dir]++;

                }
                else {
                    /* On node. Note the linear part of its (cb3, linear) bit,
                       using a reverse lookup */
                    shift_table[DECOMP_HVV_SCATTER][dir+4*index]           /* Onnode */
                        = invtab[flinear].linearcb ;
                }


                /* Gather:  mvv_recons_{plus,minus} */
                /* Operation:  chi(x) <-  \sum_dir U(x,dir)*a^F(shift(x,type=2),dir) */
                /* Receive from forward */
                if (fnode != my_node) {
                    /* Offnode */
                    /* Append to Tail 2, increase boundary count */

                    shift_table[RECONS_MVV_GATHER][dir+4*index] =
                        2*subgrid_vol_cb + (bound[cb][RECONS_MVV_GATHER][dir]);

                    bound[cb][RECONS_MVV_GATHER][dir]++;

                }
                else {
                    /* On node. Note the linear part of its (cb3, linear) bit,
                       using a reverse lookup. Note this is a recons post shift,
                       so the linear coordinate to invert is mine rather than the neighbours */
                    shift_table[RECONS_MVV_GATHER][dir+4*index] =
                        invtab[qdp_index].linearcb ;
                }

                /* Gather:  recons_{plus,minus} */
                /* Operation:  chi(x) +=  \sum_dir recons(a^B(shift(x,type=3),dir),dir) */
                /* Receive from backward */
                if (bnode != my_node) {

                    shift_table[RECONS_GATHER][dir+4*index] =
                        2*subgrid_vol_cb + bound[cb][RECONS_GATHER][dir];

                    bound[cb][RECONS_GATHER][dir]++;

                }
                else {
                    /* On node. Note the linear part of its (cb3, linear) bit,
                       using a reverse lookup. Note this is a recons post shift,
                       so the linear coordinate to invert is mine rather than the neighbours */

                    shift_table[RECONS_GATHER][dir+4*index] =
                        invtab[qdp_index].linearcb ;
                }
            }
        }
    }

    /* Sanity check - make sure the sending and receiving counters match */
    for(cb=0; cb < 2; cb++) {
        for(dir=0; dir < 4; dir++) {

            /* Sanity 1: Must have same number of boundary sites on each cb for
            a given operation */
            for(i = 0; i < 4; i++) {
                if (bound[1-cb][i][dir] != bound[cb][i][dir]) {

                    QMP_error("SSE Wilson dslash - make_shift_tables: type 0 diff. cb send/recv counts do not match: %d %d",
                              bound[1-cb][i][dir],bound[cb][i][dir]);
                    QMP_abort(1);
                }
            }


        }
    }

    /* Now I want to make the offset table into the half spinor temporaries */
    /* The half spinor temporaries will look like this:

       dir=0 [ Body Half Spinors ][ Tail 1 Half Spinors ][ Tail 2 Half Spinors ]
       dir=1 [ Body Half Spinors ][ Tail 1 Half Spinors ][ Tail 2 Half Spinors ]
       ...

       And each of these blocks of half spinors will be sized to vol_cb
       sites (ie half volume only).  The shift_table() for a given site and
       direction indexes into one of these lines. So the offset table essentially
       delineates which line one picks, by adding an offset of
       3*subgrid_vol_cb*dir
       To the shift. The result from offset table, can be used directly as a
       pointer displacement on the temporaries.

       Perhaps the best way to condsider this is to consider a value
       of shift_table[type][dir/site] that lands in the body. The
       shift table merely gives me a site index. But the data needs
       to be different for each direction for that site index. Hence
       we need to replicate the body, for each dir. The 3xsubgrid_vol_cb
       is just there to take care of the buffers.

       Or another way to think of it is that there is a 'body element' index
       specified by the shift table lookup, and that dir is just the slowest
       varying index.

    */

    /* 4 dims, 4 types, rest of the magic is to align the thingie */
    xoffset_table = (halfspinor_array **)malloc(4*4*subgrid_vol*sizeof(halfspinor_array*)+63L);
    if( xoffset_table == 0 ) {
        QMP_error("init_wnxtsu3dslash: could not initialize offset_table[i]");
        QMP_abort(1);
    }
    /* This is the bit what aligns straight from AMD Manual */
    offset_table = (halfspinor_array**)((((ptrdiff_t)(xoffset_table)) + 63L) & (-64L));

    /* Walk through the shift_table and remap the offsets into actual
       pointers */

    /* DECOMP_SCATTER */
    num=0;
    for(dir =0; dir < Nd; dir++) {

        /* Loop through all the sites. Remap the offsets either to local
           arrays or pointers */
        offsite_found=0;
        for(site=0; site < subgrid_vol; site++) {
            offset = shift_table[DECOMP_SCATTER][dir+4*site];
            if( offset >= subgrid_vol_cb ) {
                /* Found an offsite guy. It's address must be to the send back buffer */
                /* send to back index = recv from forward index = 0  */
                offsite_found++;
                offset_table[ dir + 4*(site + subgrid_vol*DECOMP_SCATTER) ] =
                    send_bufs[0][num]+(offset - subgrid_vol_cb);
            }
            else {
                /* Guy is onsite: This is DECOMP_SCATTER so offset to chi1 */
                offset_table[ dir + 4*(site + subgrid_vol*DECOMP_SCATTER) ] =
                    chi1+shift_table[DECOMP_SCATTER][dir+4*site]+subgrid_vol_cb*dir;
            }
        }

        if( offsite_found > 0 ) {
            /* If we found an offsite guy, next direction has to
            go into the next dir part of the send bufs */
            num++;
        }
    }

    /* DECOMP_HVV_SCATTER */
    /* Restart num-s */
    num=0;
    for(dir =0; dir <Nd; dir++) {
        offsite_found=0;
        for(site=0; site < subgrid_vol; site++) {
            offset = shift_table[DECOMP_HVV_SCATTER][dir+4*site];
            if( offset >= subgrid_vol_cb ) {
                /* Found an offsite guy. It's address must be to the send forw buffer */
                /* send to forward / receive from backward index = 1 */
                offsite_found++;

                offset_table[ dir + 4*(site + subgrid_vol*DECOMP_HVV_SCATTER) ] =
                    send_bufs[1][num]+(offset - subgrid_vol_cb);
            }
            else {
                /* Guy is onsite. This is DECOMP_HVV_SCATTER so offset to chi2 */
                offset_table[ dir + 4*(site + subgrid_vol*DECOMP_HVV_SCATTER) ] =
                    chi2+shift_table[DECOMP_HVV_SCATTER][dir+4*site ]+subgrid_vol_cb*dir;
            }
        }
        if( offsite_found > 0 ) {
            num++;
        }
    }

    /* RECONS_MVV_GATHER */
    num=0;
    for(dir =0; dir <Nd; dir++) {
        offsite_found=0;
        for(site=0; site < subgrid_vol; site++) {
            offset = shift_table[RECONS_MVV_GATHER][dir+4*site];
            if( offset >= 2*subgrid_vol_cb ) {
                /* Found an offsite guy. It's address must be to the recv from front buffer */
                /* recv_from front index = send to back index = 0 */
                offsite_found++;
                offset_table[ dir + 4*(site + subgrid_vol*RECONS_MVV_GATHER) ] =
                    recv_bufs[0][num]+(offset - 2*subgrid_vol_cb);
            }
            else {
                /* Guy is onsite */
                /* This is RECONS_MVV_GATHER so offset with respect to chi1 */
                offset_table[ dir + 4*(site + subgrid_vol*RECONS_MVV_GATHER) ] =
                    chi1+shift_table[RECONS_MVV_GATHER][dir+4*site ]+subgrid_vol_cb*dir;
            }
        }
        if( offsite_found > 0 ) {
            num++;
        }
    }

    /* RECONS_GATHER */
    num=0;
    for(dir =0; dir <Nd; dir++) {
        offsite_found=0;
        for(site=0; site < subgrid_vol; site++) {
            offset = shift_table[RECONS_GATHER][dir+4*site];
            if( offset >= 2*subgrid_vol_cb ) {
                /* Found an offsite guy. It's address must be to the recv from back buffer */
                /* receive from back = send to forward index =  1*/
                offsite_found++;
                offset_table[ dir + 4*(site + subgrid_vol*RECONS_GATHER) ] =
                    recv_bufs[1][num]+(offset - 2*subgrid_vol_cb);
            }
            else {
                /* Guy is onsite */
                /* This is RECONS_GATHER so offset with respect to chi2 */
                offset_table[ dir + 4*(site + subgrid_vol*RECONS_GATHER ) ] =
                    chi2+shift_table[RECONS_GATHER][dir+4*site ]+subgrid_vol_cb*dir;
            }
        }
        if( offsite_found > 0 ) {
            num++;
        }
    }



    /* Free shift table - it is no longer needed. We deal solely with offsets */
    for(i=0; i < 4; i++) {
        free( (shift_table)[i] );
    }
    free( shift_table );

    free( xinvtab );

}