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)(VecScatter ctx,Vec xin,Vec yin,InsertMode addv,ScatterMode mode)
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 & !to->use_window) {
46: /* post receives since they were not previously posted */
47: if (nrecvs) {MPI_Startall_irecv(from->starts[nrecvs]*bs,nrecvs,rwaits);}
48: }
50: #if defined(PETSC_HAVE_MPI_ALLTOALLW) && !defined(PETSC_USE_64BIT_INDICES)
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,((PetscObject)ctx)->comm);
53: } else
54: #endif
55: if (ctx->packtogether || to->use_alltoallv || to->use_window) {
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,((PetscObject)ctx)->comm);
60: #if defined(PETSC_HAVE_MPI_WIN_CREATE)
61: } else if (to->use_window) {
62: PetscInt cnt;
64: MPI_Win_fence(0,from->window);
65: for (i=0; i<nsends; i++) {
66: cnt = bs*(to->starts[i+1]-to->starts[i]);
67: MPI_Put(to->values+bs*to->starts[i],cnt,MPIU_SCALAR,to->procs[i],bs*to->winstarts[i],cnt,MPIU_SCALAR,from->window);
68: }
69: #endif
70: } else if (nsends) {
71: MPI_Startall_isend(to->starts[to->n],nsends,swaits);
72: }
73: } else {
74: /* this version packs and sends one at a time */
75: for (i=0; i<nsends; i++) {
76: PETSCMAP1(Pack)(sstarts[i+1]-sstarts[i],indices + sstarts[i],xv,svalues + bs*sstarts[i]);
77: MPI_Start_isend(sstarts[i+1]-sstarts[i],swaits+i);
78: }
79: }
81: if (!from->use_readyreceiver && to->sendfirst && !to->use_alltoallv && !to->use_window) {
82: /* post receives since they were not previously posted */
83: if (nrecvs) {MPI_Startall_irecv(from->starts[nrecvs]*bs,nrecvs,rwaits);}
84: }
85: }
87: /* take care of local scatters */
88: if (to->local.n) {
89: if (to->local.is_copy && addv == INSERT_VALUES) {
90: PetscMemcpy(yv + from->local.copy_start,xv + to->local.copy_start,to->local.copy_length);
91: } else {
92: PETSCMAP1(Scatter)(to->local.n,to->local.vslots,xv,from->local.vslots,yv,addv);
93: }
94: }
95: VecRestoreArray(xin,&xv);
96: if (xin != yin) {VecRestoreArray(yin,&yv);}
97: return(0);
98: }
100: /* --------------------------------------------------------------------------------------*/
104: PetscErrorCode PETSCMAP1(VecScatterEnd)(VecScatter ctx,Vec xin,Vec yin,InsertMode addv,ScatterMode mode)
105: {
106: VecScatter_MPI_General *to,*from;
107: PetscScalar *rvalues,*yv;
108: PetscErrorCode ierr;
109: PetscInt nrecvs,nsends,*indices,count,*rstarts,bs;
110: PetscMPIInt imdex;
111: MPI_Request *rwaits,*swaits;
112: MPI_Status xrstatus,*rstatus,*sstatus;
115: if (mode & SCATTER_LOCAL) return(0);
116: VecGetArray(yin,&yv);
118: to = (VecScatter_MPI_General*)ctx->todata;
119: from = (VecScatter_MPI_General*)ctx->fromdata;
120: rwaits = from->requests;
121: swaits = to->requests;
122: sstatus = to->sstatus; /* sstatus and rstatus are always stored in to */
123: rstatus = to->rstatus;
124: if (mode & SCATTER_REVERSE) {
125: to = (VecScatter_MPI_General*)ctx->fromdata;
126: from = (VecScatter_MPI_General*)ctx->todata;
127: rwaits = from->rev_requests;
128: swaits = to->rev_requests;
129: }
130: bs = from->bs;
131: rvalues = from->values;
132: nrecvs = from->n;
133: nsends = to->n;
134: indices = from->indices;
135: rstarts = from->starts;
137: if (ctx->packtogether || (to->use_alltoallw && (addv != INSERT_VALUES)) || (to->use_alltoallv && !to->use_alltoallw) || to->use_window) {
138: #if defined(PETSC_HAVE_MPI_WIN_CREATE)
139: if (to->use_window) {MPI_Win_fence(0,from->window);}
140: else
141: #endif
142: if (nrecvs && !to->use_alltoallv) {MPI_Waitall(nrecvs,rwaits,rstatus);}
143: PETSCMAP1(UnPack)(from->starts[from->n],from->values,indices,yv,addv);
144: } else if (!to->use_alltoallw) {
145: /* unpack one at a time */
146: count = nrecvs;
147: while (count) {
148: if (ctx->reproduce) {
149: imdex = count - 1;
150: MPI_Wait(rwaits+imdex,&xrstatus);
151: } else {
152: MPI_Waitany(nrecvs,rwaits,&imdex,&xrstatus);
153: }
154: /* unpack receives into our local space */
155: PETSCMAP1(UnPack)(rstarts[imdex+1] - rstarts[imdex],rvalues + bs*rstarts[imdex],indices + rstarts[imdex],yv,addv);
156: count--;
157: }
158: }
159: if (from->use_readyreceiver) {
160: if (nrecvs) {MPI_Startall_irecv(from->starts[nrecvs]*bs,nrecvs,rwaits);}
161: MPI_Barrier(((PetscObject)ctx)->comm);
162: }
164: /* wait on sends */
165: if (nsends && !to->use_alltoallv && !to->use_window) {MPI_Waitall(nsends,swaits,sstatus);}
166: VecRestoreArray(yin,&yv);
167: return(0);
168: }
170: #undef PETSCMAP1_a
171: #undef PETSCMAP1_b
172: #undef PETSCMAP1
173: #undef BS