Actual source code: mmaij.c

  1: /*
  2:    Support for the parallel AIJ matrix vector multiply
  3: */
  4: #include <../src/mat/impls/aij/mpi/mpiaij.h>
  5: #include <petsc/private/vecimpl.h>
  6: #include <petsc/private/isimpl.h>

  8: PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat)
  9: {
 10:   Mat_MPIAIJ *aij = (Mat_MPIAIJ *)mat->data;
 11:   Mat_SeqAIJ *B   = (Mat_SeqAIJ *)aij->B->data;
 12:   PetscInt    i, j, *aj = B->j, *garray;
 13:   PetscInt    ec = 0; /* Number of nonzero external columns */
 14:   IS          from, to;
 15:   Vec         gvec;
 16: #if defined(PETSC_USE_CTABLE)
 17:   PetscHMapI    gid1_lid1 = NULL;
 18:   PetscHashIter tpos;
 19:   PetscInt      gid, lid;
 20: #else
 21:   PetscInt N = mat->cmap->N, *indices;
 22: #endif

 24:   PetscFunctionBegin;
 25:   if (!aij->garray) {
 26:     PetscCheck(aij->B, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing B mat");
 27: #if defined(PETSC_USE_CTABLE)
 28:     /* use a table */
 29:     PetscCall(PetscHMapICreateWithSize(aij->B->rmap->n, &gid1_lid1));
 30:     for (i = 0; i < aij->B->rmap->n; i++) {
 31:       for (j = 0; j < B->ilen[i]; j++) {
 32:         PetscInt data, gid1 = aj[B->i[i] + j] + 1;
 33:         PetscCall(PetscHMapIGetWithDefault(gid1_lid1, gid1, 0, &data));
 34:         if (!data) {
 35:           /* one based table */
 36:           PetscCall(PetscHMapISet(gid1_lid1, gid1, ++ec));
 37:         }
 38:       }
 39:     }
 40:     /* form array of columns we need */
 41:     PetscCall(PetscMalloc1(ec, &garray));
 42:     PetscHashIterBegin(gid1_lid1, tpos);
 43:     while (!PetscHashIterAtEnd(gid1_lid1, tpos)) {
 44:       PetscHashIterGetKey(gid1_lid1, tpos, gid);
 45:       PetscHashIterGetVal(gid1_lid1, tpos, lid);
 46:       PetscHashIterNext(gid1_lid1, tpos);
 47:       gid--;
 48:       lid--;
 49:       garray[lid] = gid;
 50:     }
 51:     PetscCall(PetscSortInt(ec, garray)); /* sort, and rebuild */
 52:     PetscCall(PetscHMapIClear(gid1_lid1));
 53:     for (i = 0; i < ec; i++) PetscCall(PetscHMapISet(gid1_lid1, garray[i] + 1, i + 1));
 54:     /* compact out the extra columns in B */
 55:     for (i = 0; i < aij->B->rmap->n; i++) {
 56:       for (j = 0; j < B->ilen[i]; j++) {
 57:         PetscInt gid1 = aj[B->i[i] + j] + 1;
 58:         PetscCall(PetscHMapIGetWithDefault(gid1_lid1, gid1, 0, &lid));
 59:         lid--;
 60:         aj[B->i[i] + j] = lid;
 61:       }
 62:     }
 63:     PetscCall(PetscLayoutDestroy(&aij->B->cmap));
 64:     PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B), ec, ec, 1, &aij->B->cmap));
 65:     PetscCall(PetscHMapIDestroy(&gid1_lid1));
 66: #else
 67:     /* Make an array as long as the number of columns */
 68:     /* mark those columns that are in aij->B */
 69:     PetscCall(PetscCalloc1(N, &indices));
 70:     for (i = 0; i < aij->B->rmap->n; i++) {
 71:       for (j = 0; j < B->ilen[i]; j++) {
 72:         if (!indices[aj[B->i[i] + j]]) ec++;
 73:         indices[aj[B->i[i] + j]] = 1;
 74:       }
 75:     }

 77:     /* form array of columns we need */
 78:     PetscCall(PetscMalloc1(ec, &garray));
 79:     ec = 0;
 80:     for (i = 0; i < N; i++) {
 81:       if (indices[i]) garray[ec++] = i;
 82:     }

 84:     /* make indices now point into garray */
 85:     for (i = 0; i < ec; i++) indices[garray[i]] = i;

 87:     /* compact out the extra columns in B */
 88:     for (i = 0; i < aij->B->rmap->n; i++) {
 89:       for (j = 0; j < B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
 90:     }
 91:     PetscCall(PetscLayoutDestroy(&aij->B->cmap));
 92:     PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B), ec, ec, 1, &aij->B->cmap));
 93:     PetscCall(PetscFree(indices));
 94: #endif
 95:   } else {
 96:     garray = aij->garray;
 97:   }

 99:   if (!aij->lvec) {
100:     PetscCheck(aij->B, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing B mat");
101:     PetscCall(MatCreateVecs(aij->B, &aij->lvec, NULL));
102:   }
103:   PetscCall(VecGetSize(aij->lvec, &ec));

105:   /* create two temporary Index sets for build scatter gather */
106:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ec, garray, PETSC_COPY_VALUES, &from));
107:   PetscCall(ISCreateStride(PETSC_COMM_SELF, ec, 0, 1, &to));

109:   /* create temporary global vector to generate scatter context */
110:   /* This does not allocate the array's memory so is efficient */
111:   PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat), 1, mat->cmap->n, mat->cmap->N, NULL, &gvec));

113:   /* generate the scatter context */
114:   PetscCall(VecScatterDestroy(&aij->Mvctx));
115:   PetscCall(VecScatterCreate(gvec, from, aij->lvec, to, &aij->Mvctx));
116:   PetscCall(VecScatterViewFromOptions(aij->Mvctx, (PetscObject)mat, "-matmult_vecscatter_view"));
117:   aij->garray = garray;

119:   PetscCall(ISDestroy(&from));
120:   PetscCall(ISDestroy(&to));
121:   PetscCall(VecDestroy(&gvec));
122:   PetscFunctionReturn(PETSC_SUCCESS);
123: }

125: /* Disassemble the off-diag portion of the MPIAIJXxx matrix.
126:    In other words, change the B from reduced format using local col ids
127:    to expanded format using global col ids, which will make it easier to
128:    insert new nonzeros (with global col ids) into the matrix.
129:    The off-diag B determines communication in the matrix vector multiply.
130: */
131: PetscErrorCode MatDisAssemble_MPIAIJ(Mat A)
132: {
133:   Mat_MPIAIJ *aij = (Mat_MPIAIJ *)A->data;
134:   Mat         B = aij->B, Bnew = NULL;

136:   PetscFunctionBegin;
137:   /* free stuff related to matrix-vec multiply */
138:   PetscCall(VecDestroy(&aij->lvec));
139:   if (aij->colmap) {
140: #if defined(PETSC_USE_CTABLE)
141:     PetscCall(PetscHMapIDestroy(&aij->colmap));
142: #else
143:     PetscCall(PetscFree(aij->colmap));
144: #endif
145:   }

147:   if (B) {
148:     Mat_SeqAIJ        *Baij = (Mat_SeqAIJ *)B->data;
149:     PetscInt           i, j, m = B->rmap->n, n = A->cmap->N, col, ct = 0, *garray = aij->garray, *nz;
150:     PetscScalar        v;
151:     const PetscScalar *ba;

153:     /* make sure that B is assembled so we can access its values */
154:     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
155:     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));

157:     /* invent new B and copy stuff over */
158:     PetscCall(PetscMalloc1(m + 1, &nz));
159:     for (i = 0; i < m; i++) nz[i] = Baij->i[i + 1] - Baij->i[i];
160:     PetscCall(MatCreate(PETSC_COMM_SELF, &Bnew));
161:     PetscCall(MatSetSizes(Bnew, m, n, m, n)); /* Bnew now uses A->cmap->N as its col size */
162:     PetscCall(MatSetBlockSizesFromMats(Bnew, A, A));
163:     PetscCall(MatSetType(Bnew, ((PetscObject)B)->type_name));
164:     PetscCall(MatSeqAIJSetPreallocation(Bnew, 0, nz));

166:     if (Baij->nonew >= 0) { /* Inherit insertion error options (if positive). */
167:       ((Mat_SeqAIJ *)Bnew->data)->nonew = Baij->nonew;
168:     }

170:     /*
171:      Ensure that B's nonzerostate is monotonically increasing.
172:      Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
173:      */
174:     Bnew->nonzerostate = B->nonzerostate;

176:     PetscCall(PetscFree(nz));
177:     PetscCall(MatSeqAIJGetArrayRead(B, &ba));
178:     for (i = 0; i < m; i++) {
179:       for (j = Baij->i[i]; j < Baij->i[i + 1]; j++) {
180:         col = garray[Baij->j[ct]];
181:         v   = ba[ct++];
182:         PetscCall(MatSetValues(Bnew, 1, &i, 1, &col, &v, B->insertmode));
183:       }
184:     }
185:     PetscCall(MatSeqAIJRestoreArrayRead(B, &ba));
186:     PetscCall(MatDestroy(&B));
187:   }
188:   PetscCall(PetscFree(aij->garray));

190:   aij->B           = Bnew;
191:   A->was_assembled = PETSC_FALSE;
192:   A->assembled     = PETSC_FALSE;
193:   PetscFunctionReturn(PETSC_SUCCESS);
194: }

196: /*      ugly stuff added for Glenn someday we should fix this up */

198: static PetscInt *auglyrmapd = NULL, *auglyrmapo = NULL; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
199: static Vec       auglydd = NULL, auglyoo = NULL;        /* work vectors used to scale the two parts of the local matrix */

201: static PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA, Vec scale)
202: {
203:   Mat_MPIAIJ *ina = (Mat_MPIAIJ *)inA->data; /*access private part of matrix */
204:   PetscInt    i, n, nt, cstart, cend, no, *garray = ina->garray, *lindices;
205:   PetscInt   *r_rmapd, *r_rmapo;

207:   PetscFunctionBegin;
208:   PetscCall(MatGetOwnershipRange(inA, &cstart, &cend));
209:   PetscCall(MatGetSize(ina->A, NULL, &n));
210:   PetscCall(PetscCalloc1(inA->rmap->mapping->n + 1, &r_rmapd));
211:   nt = 0;
212:   for (i = 0; i < inA->rmap->mapping->n; i++) {
213:     if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
214:       nt++;
215:       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
216:     }
217:   }
218:   PetscCheck(nt == n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Hmm nt %" PetscInt_FMT " n %" PetscInt_FMT, nt, n);
219:   PetscCall(PetscMalloc1(n + 1, &auglyrmapd));
220:   for (i = 0; i < inA->rmap->mapping->n; i++) {
221:     if (r_rmapd[i]) auglyrmapd[(r_rmapd[i] - 1) - cstart] = i;
222:   }
223:   PetscCall(PetscFree(r_rmapd));
224:   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &auglydd));

226:   PetscCall(PetscCalloc1(inA->cmap->N + 1, &lindices));
227:   for (i = 0; i < ina->B->cmap->n; i++) lindices[garray[i]] = i + 1;
228:   no = inA->rmap->mapping->n - nt;
229:   PetscCall(PetscCalloc1(inA->rmap->mapping->n + 1, &r_rmapo));
230:   nt = 0;
231:   for (i = 0; i < inA->rmap->mapping->n; i++) {
232:     if (lindices[inA->rmap->mapping->indices[i]]) {
233:       nt++;
234:       r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
235:     }
236:   }
237:   PetscCheck(nt <= no, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Hmm nt %" PetscInt_FMT " no %" PetscInt_FMT, nt, n);
238:   PetscCall(PetscFree(lindices));
239:   PetscCall(PetscMalloc1(nt + 1, &auglyrmapo));
240:   for (i = 0; i < inA->rmap->mapping->n; i++) {
241:     if (r_rmapo[i]) auglyrmapo[(r_rmapo[i] - 1)] = i;
242:   }
243:   PetscCall(PetscFree(r_rmapo));
244:   PetscCall(VecCreateSeq(PETSC_COMM_SELF, nt, &auglyoo));
245:   PetscFunctionReturn(PETSC_SUCCESS);
246: }

248: PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A, Vec scale)
249: {
250:   Mat_MPIAIJ        *a = (Mat_MPIAIJ *)A->data; /*access private part of matrix */
251:   PetscInt           n, i;
252:   PetscScalar       *d, *o;
253:   const PetscScalar *s;

255:   PetscFunctionBegin;
256:   if (!auglyrmapd) PetscCall(MatMPIAIJDiagonalScaleLocalSetUp(A, scale));

258:   PetscCall(VecGetArrayRead(scale, &s));

260:   PetscCall(VecGetLocalSize(auglydd, &n));
261:   PetscCall(VecGetArray(auglydd, &d));
262:   for (i = 0; i < n; i++) { d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ }
263:   PetscCall(VecRestoreArray(auglydd, &d));
264:   /* column scale "diagonal" portion of local matrix */
265:   PetscCall(MatDiagonalScale(a->A, NULL, auglydd));

267:   PetscCall(VecGetLocalSize(auglyoo, &n));
268:   PetscCall(VecGetArray(auglyoo, &o));
269:   for (i = 0; i < n; i++) { o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ }
270:   PetscCall(VecRestoreArrayRead(scale, &s));
271:   PetscCall(VecRestoreArray(auglyoo, &o));
272:   /* column scale "off-diagonal" portion of local matrix */
273:   PetscCall(MatDiagonalScale(a->B, NULL, auglyoo));
274:   PetscFunctionReturn(PETSC_SUCCESS);
275: }