Actual source code: vpscat.h

  2: /*
  3:      Defines the methods VecScatterBegin/End_1,2,......
  4:      This is included by vpscat.c with different values for BS

  6:      This is a terrible way of doing "templates" in C.
  7: */
  8: #define PETSCMAP1_a(a,b)  a ## _ ## b
  9: #define PETSCMAP1_b(a,b)  PETSCMAP1_a(a,b)
 10: #define PETSCMAP1(a)      PETSCMAP1_b(a,BS)

 14: PetscErrorCode PETSCMAP1(VecScatterBegin)(Vec xin,Vec yin,InsertMode addv,ScatterMode mode,VecScatter ctx)
 15: {
 16:   VecScatter_MPI_General *to,*from;
 17:   PetscScalar            *xv,*yv,*svalues;
 18:   MPI_Request            *rwaits,*swaits;
 19:   PetscErrorCode         ierr;
 20:   PetscInt               i,*indices,*sstarts,nrecvs,nsends,bs;

 23:   VecGetArray(xin,&xv);
 24:   if (xin != yin) {VecGetArray(yin,&yv);} else {yv = xv;}
 25:   if (mode & SCATTER_REVERSE) {
 26:     to   = (VecScatter_MPI_General*)ctx->fromdata;
 27:     from = (VecScatter_MPI_General*)ctx->todata;
 28:     rwaits   = from->rev_requests;
 29:     swaits   = to->rev_requests;
 30:   } else {
 31:     to   = (VecScatter_MPI_General*)ctx->todata;
 32:     from = (VecScatter_MPI_General*)ctx->fromdata;
 33:     rwaits   = from->requests;
 34:     swaits   = to->requests;
 35:   }
 36:   bs       = to->bs;
 37:   svalues  = to->values;
 38:   nrecvs   = from->n;
 39:   nsends   = to->n;
 40:   indices  = to->indices;
 41:   sstarts  = to->starts;

 43:   if (!(mode & SCATTER_LOCAL)) {

 45:     if (!from->use_readyreceiver && !to->sendfirst && !to->use_alltoallv) {
 46:       /* post receives since they were not previously posted    */
 47:       if (nrecvs) {MPI_Startall_irecv(from->starts[nrecvs]*from->bs,nrecvs,rwaits);}
 48:     }

 50: #if defined(PETSC_HAVE_MPI_ALLTOALLW) 
 51:     if (to->use_alltoallw && addv == INSERT_VALUES) {
 52:       MPI_Alltoallw(xv,to->wcounts,to->wdispls,to->types,yv,from->wcounts,from->wdispls,from->types,ctx->comm);
 53:     } else
 54: #endif
 55:     if (ctx->packtogether || to->use_alltoallv) {
 56:       /* this version packs all the messages together and sends, when -vecscatter_packtogether used */
 57:       PETSCMAP1(Pack)(sstarts[nsends],indices,xv,svalues);
 58:       if (to->use_alltoallv) {
 59:         MPI_Alltoallv(to->values,to->counts,to->displs,MPIU_SCALAR,from->values,from->counts,from->displs,MPIU_SCALAR,ctx->comm);
 60:       } else if (nsends) {
 61:         MPI_Startall_isend(to->starts[to->n],nsends,swaits);
 62:       }
 63:     } else {
 64:       /* this version packs and sends one at a time */
 65:       for (i=0; i<nsends; i++) {
 66:         PETSCMAP1(Pack)(sstarts[i+1]-sstarts[i],indices + sstarts[i],xv,svalues + bs*sstarts[i]);
 67:         MPI_Start_isend(sstarts[i+1]-sstarts[i],swaits+i);
 68:       }
 69:     }

 71:     if (!from->use_readyreceiver && to->sendfirst && !to->use_alltoallv) {
 72:       /* post receives since they were not previously posted   */
 73:       if (nrecvs) {MPI_Startall_irecv(from->starts[nrecvs]*from->bs,nrecvs,rwaits);}
 74:     }
 75:   }

 77:   /* take care of local scatters */
 78:   if (to->local.n) {
 79:     if (to->local.is_copy && addv == INSERT_VALUES) {
 80:       PetscMemcpy(yv + from->local.copy_start,xv + to->local.copy_start,to->local.copy_length);
 81:     } else {
 82:       PETSCMAP1(Scatter)(to->local.n,to->local.vslots,xv,from->local.vslots,yv,addv);
 83:     }
 84:   }
 85:   VecRestoreArray(xin,&xv);
 86:   if (xin != yin) {VecRestoreArray(yin,&yv);}
 87:   return(0);
 88: }

 90: /* --------------------------------------------------------------------------------------*/

 94: PetscErrorCode PETSCMAP1(VecScatterEnd)(Vec xin,Vec yin,InsertMode addv,ScatterMode mode,VecScatter ctx)
 95: {
 96:   VecScatter_MPI_General *to,*from;
 97:   PetscScalar            *rvalues,*yv;
 98:   PetscErrorCode         ierr;
 99:   PetscInt               nrecvs,nsends,*indices,count,*rstarts,bs;
100:   PetscMPIInt            imdex;
101:   MPI_Request            *rwaits,*swaits;
102:   MPI_Status             xrstatus,*rstatus,*sstatus;

105:   if (mode & SCATTER_LOCAL) return(0);
106:   VecGetArray(yin,&yv);

108:   to       = (VecScatter_MPI_General*)ctx->todata;
109:   from     = (VecScatter_MPI_General*)ctx->fromdata;
110:   sstatus  = to->sstatus;
111:   rstatus  = to->rstatus;
112:   rwaits   = from->requests;
113:   swaits   = to->requests;
114:   if (mode & SCATTER_REVERSE) {
115:     to       = (VecScatter_MPI_General*)ctx->fromdata;
116:     from     = (VecScatter_MPI_General*)ctx->todata;
117:     rwaits   = from->rev_requests;
118:     swaits   = to->rev_requests;
119:   }
120:   bs       = from->bs;
121:   rvalues  = from->values;
122:   nrecvs   = from->n;
123:   nsends   = to->n;
124:   indices  = from->indices;
125:   rstarts  = from->starts;

127:   if (ctx->packtogether || (to->use_alltoallw && (addv != INSERT_VALUES)) || (to->use_alltoallv && !to->use_alltoallw)) {
128:     if (nrecvs && !to->use_alltoallv) {MPI_Waitall(nrecvs,rwaits,rstatus);}
129:     PETSCMAP1(UnPack)(from->starts[from->n],from->values,indices,yv,addv);
130:   } else if (!to->use_alltoallw) {
131:     /* unpack one at a time */
132:     count = nrecvs;
133:     while (count) {
134:       MPI_Waitany(nrecvs,rwaits,&imdex,&xrstatus);
135:       /* unpack receives into our local space */
136:       PETSCMAP1(UnPack)(rstarts[imdex+1] - rstarts[imdex],rvalues + bs*rstarts[imdex],indices + rstarts[imdex],yv,addv);
137:       count--;
138:     }
139:   }
140:   if (from->use_readyreceiver) {
141:     if (from->n) {MPI_Startall_irecv(from->starts[from->n]*from->bs,from->n,from->requests);}
142:     MPI_Barrier(ctx->comm);
143:   }

145:   /* wait on sends */
146:   if (nsends  && !to->use_alltoallv) {MPI_Waitall(nsends,swaits,sstatus);}
147:   VecRestoreArray(yin,&yv);
148:   return(0);
149: }

151: #undef PETSCMAP1_a
152: #undef PETSCMAP1_b
153: #undef PETSCMAP1
154: #undef BS