xref: /petsc/src/mat/impls/sell/mpi/mmsell.c (revision 5e71baeff2f3138f93cd4f5927dfd596eb8325cc)
1 /*
2    Support for the parallel SELL matrix vector multiply
3 */
4 #include <../src/mat/impls/sell/mpi/mpisell.h>
5 #include <petsc/private/isimpl.h>    /* needed because accesses data structure of ISLocalToGlobalMapping directly */
6 
7 /*
8    Takes the local part of an already assembled MPISELL matrix
9    and disassembles it. This is to allow new nonzeros into the matrix
10    that require more communication in the matrix vector multiply.
11    Thus certain data-structures must be rebuilt.
12 
13    Kind of slow! But that's what application programmers get when
14    they are sloppy.
15 */
16 PetscErrorCode MatDisAssemble_MPISELL(Mat A)
17 {
18   Mat_MPISELL    *sell=(Mat_MPISELL*)A->data;
19   Mat            B=sell->B,Bnew;
20   Mat_SeqSELL    *Bsell=(Mat_SeqSELL*)B->data;
21   PetscInt       i,j,totalslices,N=A->cmap->N,ec,row;
22   PetscBool      isnonzero;
23   PetscErrorCode ierr;
24 
25   PetscFunctionBegin;
26   /* free stuff related to matrix-vec multiply */
27   ierr = VecGetSize(sell->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
28   ierr = VecDestroy(&sell->lvec);CHKERRQ(ierr);
29   ierr = VecScatterDestroy(&sell->Mvctx);CHKERRQ(ierr);
30   if (sell->colmap) {
31 #if defined(PETSC_USE_CTABLE)
32     ierr = PetscTableDestroy(&sell->colmap);CHKERRQ(ierr);
33 #else
34     ierr = PetscFree(sell->colmap);CHKERRQ(ierr);
35     ierr = PetscLogObjectMemory((PetscObject)A,-sell->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr);
36 #endif
37   }
38 
39   /* make sure that B is assembled so we can access its values */
40   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
41   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
42 
43   /* invent new B and copy stuff over */
44   ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr);
45   ierr = MatSetSizes(Bnew,B->rmap->n,N,B->rmap->n,N);CHKERRQ(ierr);
46   ierr = MatSetBlockSizesFromMats(Bnew,A,A);CHKERRQ(ierr);
47   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
48   ierr = MatSeqSELLSetPreallocation(Bnew,0,Bsell->rlen);CHKERRQ(ierr);
49   ((Mat_SeqSELL*)Bnew->data)->nonew = Bsell->nonew; /* Inherit insertion error options. */
50 
51   /*
52    Ensure that B's nonzerostate is monotonically increasing.
53    Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
54    */
55   Bnew->nonzerostate = B->nonzerostate;
56 
57   totalslices = B->rmap->n/8+((B->rmap->n & 0x07)?1:0); /* floor(n/8) */
58   for (i=0; i<totalslices; i++) { /* loop over slices */
59     for (j=Bsell->sliidx[i],row=0; j<Bsell->sliidx[i+1]; j++,row=((row+1)&0x07)) {
60       isnonzero = (PetscBool)((j-Bsell->sliidx[i])/8 < Bsell->rlen[8*i+row]);
61       if (isnonzero) {
62         ierr = MatSetValue(Bnew,8*i+row,sell->garray[Bsell->colidx[j]],Bsell->val[j],B->insertmode);CHKERRQ(ierr);
63       }
64     }
65   }
66 
67   ierr = PetscFree(sell->garray);CHKERRQ(ierr);
68   ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
69   ierr = MatDestroy(&B);CHKERRQ(ierr);
70   ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr);
71 
72   sell->B          = Bnew;
73   A->was_assembled = PETSC_FALSE;
74   A->assembled     = PETSC_FALSE;
75   PetscFunctionReturn(0);
76 }
77 
78 PetscErrorCode MatSetUpMultiply_MPISELL(Mat mat)
79 {
80   Mat_MPISELL    *sell=(Mat_MPISELL*)mat->data;
81   Mat_SeqSELL    *B=(Mat_SeqSELL*)(sell->B->data);
82   PetscErrorCode ierr;
83   PetscInt       i,j,*bcolidx=B->colidx,ec=0,*garray,totalslices;
84   IS             from,to;
85   Vec            gvec;
86   PetscBool      isnonzero;
87 #if defined(PETSC_USE_CTABLE)
88   PetscTable         gid1_lid1;
89   PetscTablePosition tpos;
90   PetscInt           gid,lid;
91 #else
92   PetscInt N = mat->cmap->N,*indices;
93 #endif
94 
95   PetscFunctionBegin;
96   totalslices = sell->B->rmap->n/8+((sell->B->rmap->n & 0x07)?1:0); /* floor(n/8) */
97 
98   /* ec counts the number of columns that contain nonzeros */
99 #if defined(PETSC_USE_CTABLE)
100   /* use a table */
101   ierr = PetscTableCreate(sell->B->rmap->n,mat->cmap->N+1,&gid1_lid1);CHKERRQ(ierr);
102   for (i=0; i<totalslices; i++) { /* loop over slices */
103     for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) {
104       isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]);
105       if (isnonzero) { /* check the mask bit */
106         PetscInt data,gid1 = bcolidx[j] + 1;
107         ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
108         if (!data) {
109           /* one based table */
110           ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr);
111         }
112       }
113     }
114   }
115 
116   /* form array of columns we need */
117   ierr = PetscMalloc1(ec+1,&garray);CHKERRQ(ierr);
118   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
119   while (tpos) {
120     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
121     gid--;
122     lid--;
123     garray[lid] = gid;
124   }
125   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */
126   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
127   for (i=0; i<ec; i++) {
128     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr);
129   }
130 
131   /* compact out the extra columns in B */
132   for (i=0; i<totalslices; i++) { /* loop over slices */
133     for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) {
134       isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]);
135       if (isnonzero) {
136         PetscInt gid1 = bcolidx[j] + 1;
137         ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
138         lid--;
139         bcolidx[j] = lid;
140       }
141     }
142   }
143   sell->B->cmap->n = sell->B->cmap->N = ec;
144   sell->B->cmap->bs = 1;
145 
146   ierr = PetscLayoutSetUp((sell->B->cmap));CHKERRQ(ierr);
147   ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr);
148 #else
149   /* Make an array as long as the number of columns */
150   ierr = PetscCalloc1(N+1,&indices);CHKERRQ(ierr);
151   /* mark those columns that are in sell->B */
152   for (i=0; i<totalslices; i++) { /* loop over slices */
153     for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) {
154       isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]);
155       if (isnonzero) {
156         if (!indices[bcolidx[j]]) ec++;
157         indices[bcolidx[j]] = 1;
158       }
159     }
160   }
161 
162   /* form array of columns we need */
163   ierr = PetscMalloc1(ec+1,&garray);CHKERRQ(ierr);
164   ec   = 0;
165   for (i=0; i<N; i++) {
166     if (indices[i]) garray[ec++] = i;
167   }
168 
169   /* make indices now point into garray */
170   for (i=0; i<ec; i++) {
171     indices[garray[i]] = i;
172   }
173 
174   /* compact out the extra columns in B */
175   for (i=0; i<totalslices; i++) { /* loop over slices */
176     for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) {
177       isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]);
178       if (isnonzero) bcolidx[j] = indices[bcolidx[j]];
179     }
180   }
181   sell->B->cmap->n = sell->B->cmap->N = ec; /* number of columns that are not all zeros */
182   sell->B->cmap->bs = 1;
183 
184   ierr = PetscLayoutSetUp((sell->B->cmap));CHKERRQ(ierr);
185   ierr = PetscFree(indices);CHKERRQ(ierr);
186 #endif
187   /* create local vector that is used to scatter into */
188   ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&sell->lvec);CHKERRQ(ierr);
189   /* create two temporary Index sets for build scatter gather */
190   ierr = ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
191   ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr);
192 
193   /* create temporary global vector to generate scatter context */
194   /* This does not allocate the array's memory so is efficient */
195   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
196 
197   /* generate the scatter context */
198   ierr = VecScatterCreate(gvec,from,sell->lvec,to,&sell->Mvctx);CHKERRQ(ierr);
199   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sell->Mvctx);CHKERRQ(ierr);
200   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sell->lvec);CHKERRQ(ierr);
201   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
202   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
203 
204   sell->garray = garray;
205 
206   ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
207   ierr = ISDestroy(&from);CHKERRQ(ierr);
208   ierr = ISDestroy(&to);CHKERRQ(ierr);
209   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
210   PetscFunctionReturn(0);
211 }
212 
213 /*      ugly stuff added for Glenn someday we should fix this up */
214 static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
215 static Vec auglydd          = 0,auglyoo     = 0; /* work vectors used to scale the two parts of the local matrix */
216 
217 PetscErrorCode MatMPISELLDiagonalScaleLocalSetUp(Mat inA,Vec scale)
218 {
219   Mat_MPISELL    *ina=(Mat_MPISELL*)inA->data; /*access private part of matrix */
220   PetscErrorCode ierr;
221   PetscInt       i,n,nt,cstart,cend,no,*garray=ina->garray,*lindices;
222   PetscInt       *r_rmapd,*r_rmapo;
223 
224   PetscFunctionBegin;
225   ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr);
226   ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr);
227   ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr);
228   nt   = 0;
229   for (i=0; i<inA->rmap->mapping->n; i++) {
230     if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
231       nt++;
232       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
233     }
234   }
235   if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n);
236   ierr = PetscMalloc1(n+1,&auglyrmapd);CHKERRQ(ierr);
237   for (i=0; i<inA->rmap->mapping->n; i++) {
238     if (r_rmapd[i]) {
239       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
240     }
241   }
242   ierr = PetscFree(r_rmapd);CHKERRQ(ierr);
243   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr);
244   ierr = PetscCalloc1(inA->cmap->N+1,&lindices);CHKERRQ(ierr);
245   for (i=0; i<ina->B->cmap->n; i++) {
246     lindices[garray[i]] = i+1;
247   }
248   no   = inA->rmap->mapping->n - nt;
249   ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr);
250   nt   = 0;
251   for (i=0; i<inA->rmap->mapping->n; i++) {
252     if (lindices[inA->rmap->mapping->indices[i]]) {
253       nt++;
254       r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
255     }
256   }
257   if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
258   ierr = PetscFree(lindices);CHKERRQ(ierr);
259   ierr = PetscMalloc1(nt+1,&auglyrmapo);CHKERRQ(ierr);
260   for (i=0; i<inA->rmap->mapping->n; i++) {
261     if (r_rmapo[i]) {
262       auglyrmapo[(r_rmapo[i]-1)] = i;
263     }
264   }
265   ierr = PetscFree(r_rmapo);CHKERRQ(ierr);
266   ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr);
267   PetscFunctionReturn(0);
268 }
269 
270 PetscErrorCode MatDiagonalScaleLocal_MPISELL(Mat A,Vec scale)
271 {
272   Mat_MPISELL       *a=(Mat_MPISELL*)A->data; /*access private part of matrix */
273   PetscErrorCode    ierr;
274   PetscInt          n,i;
275   PetscScalar       *d,*o;
276   const PetscScalar *s;
277 
278   PetscFunctionBegin;
279   if (!auglyrmapd) {
280     ierr = MatMPISELLDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr);
281   }
282   ierr = VecGetArrayRead(scale,&s);CHKERRQ(ierr);
283   ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr);
284   ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr);
285   for (i=0; i<n; i++) {
286     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
287   }
288   ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr);
289   /* column scale "diagonal" portion of local matrix */
290   ierr = MatDiagonalScale(a->A,NULL,auglydd);CHKERRQ(ierr);
291   ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr);
292   ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr);
293   for (i=0; i<n; i++) {
294     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
295   }
296   ierr = VecRestoreArrayRead(scale,&s);CHKERRQ(ierr);
297   ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr);
298   /* column scale "off-diagonal" portion of local matrix */
299   ierr = MatDiagonalScale(a->B,NULL,auglyoo);CHKERRQ(ierr);
300   PetscFunctionReturn(0);
301 }
302