xref: /petsc/src/mat/impls/aij/mpi/mmaij.c (revision b8d709ab44c354dba7591c0191c0984ef3bb0cd5)
1 
2 /*
3    Support for the parallel AIJ matrix vector multiply
4 */
5 #include <../src/mat/impls/aij/mpi/mpiaij.h>
6 #include <petsc/private/vecimpl.h>
7 #include <petsc/private/isimpl.h> /* needed because accesses data structure of ISLocalToGlobalMapping directly */
8 
9 PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat)
10 {
11   Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data;
12   Mat_SeqAIJ *B   = (Mat_SeqAIJ *)(aij->B->data);
13   PetscInt    i, j, *aj = B->j, *garray;
14   PetscInt    ec = 0; /* Number of nonzero external columns */
15   IS          from, to;
16   Vec         gvec;
17 #if defined(PETSC_USE_CTABLE)
18   PetscTable         gid1_lid1;
19   PetscTablePosition tpos;
20   PetscInt           gid, lid;
21 #else
22   PetscInt N = mat->cmap->N, *indices;
23 #endif
24 
25   PetscFunctionBegin;
26   if (!aij->garray) {
27     PetscCheck(aij->B, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing B mat");
28 #if defined(PETSC_USE_CTABLE)
29     /* use a table */
30     PetscCall(PetscTableCreate(aij->B->rmap->n, mat->cmap->N + 1, &gid1_lid1));
31     for (i = 0; i < aij->B->rmap->n; i++) {
32       for (j = 0; j < B->ilen[i]; j++) {
33         PetscInt data, gid1 = aj[B->i[i] + j] + 1;
34         PetscCall(PetscTableFind(gid1_lid1, gid1, &data));
35         if (!data) {
36           /* one based table */
37           PetscCall(PetscTableAdd(gid1_lid1, gid1, ++ec, INSERT_VALUES));
38         }
39       }
40     }
41     /* form array of columns we need */
42     PetscCall(PetscMalloc1(ec, &garray));
43     PetscCall(PetscTableGetHeadPosition(gid1_lid1, &tpos));
44     while (tpos) {
45       PetscCall(PetscTableGetNext(gid1_lid1, &tpos, &gid, &lid));
46       gid--;
47       lid--;
48       garray[lid] = gid;
49     }
50     PetscCall(PetscSortInt(ec, garray)); /* sort, and rebuild */
51     PetscCall(PetscTableRemoveAll(gid1_lid1));
52     for (i = 0; i < ec; i++) PetscCall(PetscTableAdd(gid1_lid1, garray[i] + 1, i + 1, INSERT_VALUES));
53     /* compact out the extra columns in B */
54     for (i = 0; i < aij->B->rmap->n; i++) {
55       for (j = 0; j < B->ilen[i]; j++) {
56         PetscInt gid1 = aj[B->i[i] + j] + 1;
57         PetscCall(PetscTableFind(gid1_lid1, gid1, &lid));
58         lid--;
59         aj[B->i[i] + j] = lid;
60       }
61     }
62     PetscCall(PetscLayoutDestroy(&aij->B->cmap));
63     PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B), ec, ec, 1, &aij->B->cmap));
64     PetscCall(PetscTableDestroy(&gid1_lid1));
65 #else
66     /* Make an array as long as the number of columns */
67     /* mark those columns that are in aij->B */
68     PetscCall(PetscCalloc1(N, &indices));
69     for (i = 0; i < aij->B->rmap->n; i++) {
70       for (j = 0; j < B->ilen[i]; j++) {
71         if (!indices[aj[B->i[i] + j]]) ec++;
72         indices[aj[B->i[i] + j]] = 1;
73       }
74     }
75 
76     /* form array of columns we need */
77     PetscCall(PetscMalloc1(ec, &garray));
78     ec = 0;
79     for (i = 0; i < N; i++) {
80       if (indices[i]) garray[ec++] = i;
81     }
82 
83     /* make indices now point into garray */
84     for (i = 0; i < ec; i++) indices[garray[i]] = i;
85 
86     /* compact out the extra columns in B */
87     for (i = 0; i < aij->B->rmap->n; i++) {
88       for (j = 0; j < B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
89     }
90     PetscCall(PetscLayoutDestroy(&aij->B->cmap));
91     PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B), ec, ec, 1, &aij->B->cmap));
92     PetscCall(PetscFree(indices));
93 #endif
94   } else {
95     garray = aij->garray;
96   }
97 
98   if (!aij->lvec) {
99     PetscCheck(aij->B, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing B mat");
100     PetscCall(MatCreateVecs(aij->B, &aij->lvec, NULL));
101   }
102   PetscCall(VecGetSize(aij->lvec, &ec));
103 
104   /* create two temporary Index sets for build scatter gather */
105   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ec, garray, PETSC_COPY_VALUES, &from));
106   PetscCall(ISCreateStride(PETSC_COMM_SELF, ec, 0, 1, &to));
107 
108   /* create temporary global vector to generate scatter context */
109   /* This does not allocate the array's memory so is efficient */
110   PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat), 1, mat->cmap->n, mat->cmap->N, NULL, &gvec));
111 
112   /* generate the scatter context */
113   PetscCall(VecScatterDestroy(&aij->Mvctx));
114   PetscCall(VecScatterCreate(gvec, from, aij->lvec, to, &aij->Mvctx));
115   PetscCall(VecScatterViewFromOptions(aij->Mvctx, (PetscObject)mat, "-matmult_vecscatter_view"));
116   aij->garray = garray;
117 
118   PetscCall(ISDestroy(&from));
119   PetscCall(ISDestroy(&to));
120   PetscCall(VecDestroy(&gvec));
121   PetscFunctionReturn(0);
122 }
123 
124 /* Disassemble the off-diag portion of the MPIAIJXxx matrix.
125    In other words, change the B from reduced format using local col ids
126    to expanded format using global col ids, which will make it easier to
127    insert new nonzeros (with global col ids) into the matrix.
128    The off-diag B determines communication in the matrix vector multiply.
129 */
130 PetscErrorCode MatDisAssemble_MPIAIJ(Mat A)
131 {
132   Mat_MPIAIJ        *aij  = (Mat_MPIAIJ *)A->data;
133   Mat                B    = aij->B, Bnew;
134   Mat_SeqAIJ        *Baij = (Mat_SeqAIJ *)B->data;
135   PetscInt           i, j, m = B->rmap->n, n = A->cmap->N, col, ct = 0, *garray = aij->garray, *nz;
136   PetscScalar        v;
137   const PetscScalar *ba;
138 
139   PetscFunctionBegin;
140   /* free stuff related to matrix-vec multiply */
141   PetscCall(VecDestroy(&aij->lvec));
142   if (aij->colmap) {
143 #if defined(PETSC_USE_CTABLE)
144     PetscCall(PetscTableDestroy(&aij->colmap));
145 #else
146     PetscCall(PetscFree(aij->colmap));
147 #endif
148   }
149 
150   /* make sure that B is assembled so we can access its values */
151   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
152   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
153 
154   /* invent new B and copy stuff over */
155   PetscCall(PetscMalloc1(m + 1, &nz));
156   for (i = 0; i < m; i++) nz[i] = Baij->i[i + 1] - Baij->i[i];
157   PetscCall(MatCreate(PETSC_COMM_SELF, &Bnew));
158   PetscCall(MatSetSizes(Bnew, m, n, m, n)); /* Bnew now uses A->cmap->N as its col size */
159   PetscCall(MatSetBlockSizesFromMats(Bnew, A, A));
160   PetscCall(MatSetType(Bnew, ((PetscObject)B)->type_name));
161   PetscCall(MatSeqAIJSetPreallocation(Bnew, 0, nz));
162 
163   if (Baij->nonew >= 0) { /* Inherit insertion error options (if positive). */
164     ((Mat_SeqAIJ *)Bnew->data)->nonew = Baij->nonew;
165   }
166 
167   /*
168    Ensure that B's nonzerostate is monotonically increasing.
169    Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
170    */
171   Bnew->nonzerostate = B->nonzerostate;
172 
173   PetscCall(PetscFree(nz));
174   PetscCall(MatSeqAIJGetArrayRead(B, &ba));
175   for (i = 0; i < m; i++) {
176     for (j = Baij->i[i]; j < Baij->i[i + 1]; j++) {
177       col = garray[Baij->j[ct]];
178       v   = ba[ct++];
179       PetscCall(MatSetValues(Bnew, 1, &i, 1, &col, &v, B->insertmode));
180     }
181   }
182   PetscCall(MatSeqAIJRestoreArrayRead(B, &ba));
183 
184   PetscCall(PetscFree(aij->garray));
185   PetscCall(MatDestroy(&B));
186 
187   aij->B           = Bnew;
188   A->was_assembled = PETSC_FALSE;
189   A->assembled     = PETSC_FALSE;
190   PetscFunctionReturn(0);
191 }
192 
193 /*      ugly stuff added for Glenn someday we should fix this up */
194 
195 static PetscInt *auglyrmapd = NULL, *auglyrmapo = NULL; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
196 static Vec       auglydd = NULL, auglyoo = NULL;        /* work vectors used to scale the two parts of the local matrix */
197 
198 PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA, Vec scale)
199 {
200   Mat_MPIAIJ *ina = (Mat_MPIAIJ *)inA->data; /*access private part of matrix */
201   PetscInt    i, n, nt, cstart, cend, no, *garray = ina->garray, *lindices;
202   PetscInt   *r_rmapd, *r_rmapo;
203 
204   PetscFunctionBegin;
205   PetscCall(MatGetOwnershipRange(inA, &cstart, &cend));
206   PetscCall(MatGetSize(ina->A, NULL, &n));
207   PetscCall(PetscCalloc1(inA->rmap->mapping->n + 1, &r_rmapd));
208   nt = 0;
209   for (i = 0; i < inA->rmap->mapping->n; i++) {
210     if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
211       nt++;
212       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
213     }
214   }
215   PetscCheck(nt == n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Hmm nt %" PetscInt_FMT " n %" PetscInt_FMT, nt, n);
216   PetscCall(PetscMalloc1(n + 1, &auglyrmapd));
217   for (i = 0; i < inA->rmap->mapping->n; i++) {
218     if (r_rmapd[i]) auglyrmapd[(r_rmapd[i] - 1) - cstart] = i;
219   }
220   PetscCall(PetscFree(r_rmapd));
221   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &auglydd));
222 
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 MatMPIAIJDiagonalScaleLocal(Mat A, Vec scale)
246 {
247   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
248 
249   PetscFunctionBegin;
250   PetscTryMethod(A, "MatDiagonalScaleLocal_C", (Mat, Vec), (A, scale));
251   PetscFunctionReturn(0);
252 }
253 
254 PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A, Vec scale)
255 {
256   Mat_MPIAIJ        *a = (Mat_MPIAIJ *)A->data; /*access private part of matrix */
257   PetscInt           n, i;
258   PetscScalar       *d, *o;
259   const PetscScalar *s;
260 
261   PetscFunctionBegin;
262   if (!auglyrmapd) PetscCall(MatMPIAIJDiagonalScaleLocalSetUp(A, scale));
263 
264   PetscCall(VecGetArrayRead(scale, &s));
265 
266   PetscCall(VecGetLocalSize(auglydd, &n));
267   PetscCall(VecGetArray(auglydd, &d));
268   for (i = 0; i < n; i++) { d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ }
269   PetscCall(VecRestoreArray(auglydd, &d));
270   /* column scale "diagonal" portion of local matrix */
271   PetscCall(MatDiagonalScale(a->A, NULL, auglydd));
272 
273   PetscCall(VecGetLocalSize(auglyoo, &n));
274   PetscCall(VecGetArray(auglyoo, &o));
275   for (i = 0; i < n; i++) { o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ }
276   PetscCall(VecRestoreArrayRead(scale, &s));
277   PetscCall(VecRestoreArray(auglyoo, &o));
278   /* column scale "off-diagonal" portion of local matrix */
279   PetscCall(MatDiagonalScale(a->B, NULL, auglyoo));
280   PetscFunctionReturn(0);
281 }
282