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