Actual source code: vpscat.h
petsc-3.6.4 2016-04-12
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: if (mode & SCATTER_REVERSE) {
24: to = (VecScatter_MPI_General*)ctx->fromdata;
25: from = (VecScatter_MPI_General*)ctx->todata;
26: rwaits = from->rev_requests;
27: swaits = to->rev_requests;
28: } else {
29: to = (VecScatter_MPI_General*)ctx->todata;
30: from = (VecScatter_MPI_General*)ctx->fromdata;
31: rwaits = from->requests;
32: swaits = to->requests;
33: }
34: bs = to->bs;
35: svalues = to->values;
36: nrecvs = from->n;
37: nsends = to->n;
38: indices = to->indices;
39: sstarts = to->starts;
40: #if defined(PETSC_HAVE_CUSP)
41: VecCUSPAllocateCheckHost(xin);
42: if (xin->valid_GPU_array == PETSC_CUSP_GPU) {
43: if (xin->spptr && ctx->spptr) {
44: VecCUSPCopyFromGPUSome_Public(xin,(PetscCUSPIndices)ctx->spptr);
45: } else {
46: VecCUSPCopyFromGPU(xin);
47: }
48: }
49: xv = *((PetscScalar**)xin->data);
50: #else
51: VecGetArrayRead(xin,(const PetscScalar**)&xv);
52: #endif
54: if (xin != yin) {VecGetArray(yin,&yv);}
55: else yv = xv;
57: if (!(mode & SCATTER_LOCAL)) {
58: if (!from->use_readyreceiver && !to->sendfirst && !to->use_alltoallv & !to->use_window) {
59: /* post receives since they were not previously posted */
60: if (nrecvs) {MPI_Startall_irecv(from->starts[nrecvs]*bs,nrecvs,rwaits);}
61: }
63: #if defined(PETSC_HAVE_MPI_ALLTOALLW) && !defined(PETSC_USE_64BIT_INDICES)
64: if (to->use_alltoallw && addv == INSERT_VALUES) {
65: MPI_Alltoallw(xv,to->wcounts,to->wdispls,to->types,yv,from->wcounts,from->wdispls,from->types,PetscObjectComm((PetscObject)ctx));
66: } else
67: #endif
68: if (ctx->packtogether || to->use_alltoallv || to->use_window) {
69: /* this version packs all the messages together and sends, when -vecscatter_packtogether used */
70: PETSCMAP1(Pack)(sstarts[nsends],indices,xv,svalues,bs);
71: if (to->use_alltoallv) {
72: MPI_Alltoallv(to->values,to->counts,to->displs,MPIU_SCALAR,from->values,from->counts,from->displs,MPIU_SCALAR,PetscObjectComm((PetscObject)ctx));
73: #if defined(PETSC_HAVE_MPI_WIN_CREATE)
74: } else if (to->use_window) {
75: PetscInt cnt;
77: MPI_Win_fence(0,from->window);
78: for (i=0; i<nsends; i++) {
79: cnt = bs*(to->starts[i+1]-to->starts[i]);
80: MPI_Put(to->values+bs*to->starts[i],cnt,MPIU_SCALAR,to->procs[i],bs*to->winstarts[i],cnt,MPIU_SCALAR,from->window);
81: }
82: #endif
83: } else if (nsends) {
84: MPI_Startall_isend(to->starts[to->n],nsends,swaits);
85: }
86: } else {
87: /* this version packs and sends one at a time */
88: for (i=0; i<nsends; i++) {
89: PETSCMAP1(Pack)(sstarts[i+1]-sstarts[i],indices + sstarts[i],xv,svalues + bs*sstarts[i],bs);
90: MPI_Start_isend(sstarts[i+1]-sstarts[i],swaits+i);
91: }
92: }
94: if (!from->use_readyreceiver && to->sendfirst && !to->use_alltoallv && !to->use_window) {
95: /* post receives since they were not previously posted */
96: if (nrecvs) {MPI_Startall_irecv(from->starts[nrecvs]*bs,nrecvs,rwaits);}
97: }
98: }
100: /* take care of local scatters */
101: if (to->local.n) {
102: if (to->local.is_copy && addv == INSERT_VALUES) {
103: if (yv != xv || from->local.copy_start != to->local.copy_start) {
104: PetscMemcpy(yv + from->local.copy_start,xv + to->local.copy_start,to->local.copy_length);
105: }
106: } else {
107: if (xv == yv && addv == INSERT_VALUES && to->local.nonmatching_computed) {
108: /* only copy entries that do not share identical memory locations */
109: PETSCMAP1(Scatter)(to->local.n_nonmatching,to->local.slots_nonmatching,xv,from->local.slots_nonmatching,yv,addv,bs);
110: } else {
111: PETSCMAP1(Scatter)(to->local.n,to->local.vslots,xv,from->local.vslots,yv,addv,bs);
112: }
113: }
114: }
115: #if defined(PETSC_HAVE_CUSP)
116: if (xin->valid_GPU_array != PETSC_CUSP_GPU) {
117: VecRestoreArrayRead(xin,(const PetscScalar**)&xv);
118: } else {
119: VecRestoreArrayRead(xin,(const PetscScalar**)&xv);
120: }
121: #else
122: VecRestoreArrayRead(xin,(const PetscScalar**)&xv);
123: #endif
124: if (xin != yin) {VecRestoreArray(yin,&yv);}
125: return(0);
126: }
128: /* --------------------------------------------------------------------------------------*/
132: PetscErrorCode PETSCMAP1(VecScatterEnd)(VecScatter ctx,Vec xin,Vec yin,InsertMode addv,ScatterMode mode)
133: {
134: VecScatter_MPI_General *to,*from;
135: PetscScalar *rvalues,*yv;
136: PetscErrorCode ierr;
137: PetscInt nrecvs,nsends,*indices,count,*rstarts,bs;
138: PetscMPIInt imdex;
139: MPI_Request *rwaits,*swaits;
140: MPI_Status xrstatus,*rstatus,*sstatus;
143: if (mode & SCATTER_LOCAL) return(0);
144: VecGetArray(yin,&yv);
146: to = (VecScatter_MPI_General*)ctx->todata;
147: from = (VecScatter_MPI_General*)ctx->fromdata;
148: rwaits = from->requests;
149: swaits = to->requests;
150: sstatus = to->sstatus; /* sstatus and rstatus are always stored in to */
151: rstatus = to->rstatus;
152: if (mode & SCATTER_REVERSE) {
153: to = (VecScatter_MPI_General*)ctx->fromdata;
154: from = (VecScatter_MPI_General*)ctx->todata;
155: rwaits = from->rev_requests;
156: swaits = to->rev_requests;
157: }
158: bs = from->bs;
159: rvalues = from->values;
160: nrecvs = from->n;
161: nsends = to->n;
162: indices = from->indices;
163: rstarts = from->starts;
165: if (ctx->packtogether || (to->use_alltoallw && (addv != INSERT_VALUES)) || (to->use_alltoallv && !to->use_alltoallw) || to->use_window) {
166: #if defined(PETSC_HAVE_MPI_WIN_CREATE)
167: if (to->use_window) {MPI_Win_fence(0,from->window);}
168: else
169: #endif
170: if (nrecvs && !to->use_alltoallv) {MPI_Waitall(nrecvs,rwaits,rstatus);}
171: PETSCMAP1(UnPack)(from->starts[from->n],from->values,indices,yv,addv,bs);
172: } else if (!to->use_alltoallw) {
173: /* unpack one at a time */
174: count = nrecvs;
175: while (count) {
176: if (ctx->reproduce) {
177: imdex = count - 1;
178: MPI_Wait(rwaits+imdex,&xrstatus);
179: } else {
180: MPI_Waitany(nrecvs,rwaits,&imdex,&xrstatus);
181: }
182: /* unpack receives into our local space */
183: PETSCMAP1(UnPack)(rstarts[imdex+1] - rstarts[imdex],rvalues + bs*rstarts[imdex],indices + rstarts[imdex],yv,addv,bs);
184: count--;
185: }
186: }
187: if (from->use_readyreceiver) {
188: if (nrecvs) {MPI_Startall_irecv(from->starts[nrecvs]*bs,nrecvs,rwaits);}
189: MPI_Barrier(PetscObjectComm((PetscObject)ctx));
190: }
192: /* wait on sends */
193: if (nsends && !to->use_alltoallv && !to->use_window) {MPI_Waitall(nsends,swaits,sstatus);}
194: VecRestoreArray(yin,&yv);
195: return(0);
196: }
198: #undef PETSCMAP1_a
199: #undef PETSCMAP1_b
200: #undef PETSCMAP1
201: #undef BS