1 2 /* 3 Routines to compute overlapping regions of a parallel MPI matrix 4 and to find submatrices that were shared across processors. 5 */ 6 #include <../src/mat/impls/baij/mpi/mpibaij.h> 7 #include <petscbt.h> 8 9 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat, PetscInt, char **, PetscInt *, PetscInt **); 10 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat, PetscInt, PetscInt **, PetscInt **, PetscInt *); 11 extern PetscErrorCode MatGetRow_MPIBAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 12 extern PetscErrorCode MatRestoreRow_MPIBAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 13 14 PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C, PetscInt imax, IS is[], PetscInt ov) 15 { 16 PetscInt i, N = C->cmap->N, bs = C->rmap->bs, n; 17 const PetscInt *idx; 18 IS *is_new; 19 20 PetscFunctionBegin; 21 PetscCall(PetscMalloc1(imax, &is_new)); 22 /* Convert the indices into block format */ 23 PetscCall(ISCompressIndicesGeneral(N, C->rmap->n, bs, imax, is, is_new)); 24 PetscCheck(ov >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap specified"); 25 for (i = 0; i < ov; ++i) PetscCall(MatIncreaseOverlap_MPIBAIJ_Once(C, imax, is_new)); 26 for (i = 0; i < imax; i++) { 27 PetscCall(ISDestroy(&is[i])); 28 PetscCall(ISGetLocalSize(is_new[i], &n)); 29 PetscCall(ISGetIndices(is_new[i], &idx)); 30 PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)is_new[i]), bs, n, idx, PETSC_COPY_VALUES, &is[i])); 31 PetscCall(ISDestroy(&is_new[i])); 32 } 33 PetscCall(PetscFree(is_new)); 34 PetscFunctionReturn(PETSC_SUCCESS); 35 } 36 37 /* 38 Sample message format: 39 If a processor A wants processor B to process some elements corresponding 40 to index sets is[1], is[5] 41 mesg [0] = 2 (no of index sets in the mesg) 42 ----------- 43 mesg [1] = 1 => is[1] 44 mesg [2] = sizeof(is[1]); 45 ----------- 46 mesg [5] = 5 => is[5] 47 mesg [6] = sizeof(is[5]); 48 ----------- 49 mesg [7] 50 mesg [n] data(is[1]) 51 ----------- 52 mesg[n+1] 53 mesg[m] data(is[5]) 54 ----------- 55 56 Notes: 57 nrqs - no of requests sent (or to be sent out) 58 nrqr - no of requests received (which have to be or which have been processed) 59 */ 60 PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C, PetscInt imax, IS is[]) 61 { 62 Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 63 const PetscInt **idx, *idx_i; 64 PetscInt *n, *w3, *w4, **data, len; 65 PetscMPIInt size, rank, tag1, tag2, *w2, *w1, nrqr; 66 PetscInt Mbs, i, j, k, **rbuf, row, nrqs, msz, **outdat, **ptr; 67 PetscInt *ctr, *pa, *tmp, *isz, *isz1, **xdata, **rbuf2, *d_p; 68 PetscMPIInt *onodes1, *olengths1, *onodes2, *olengths2, proc = -1; 69 PetscBT *table; 70 MPI_Comm comm, *iscomms; 71 MPI_Request *s_waits1, *r_waits1, *s_waits2, *r_waits2; 72 char *t_p; 73 74 PetscFunctionBegin; 75 PetscCall(PetscObjectGetComm((PetscObject)C, &comm)); 76 size = c->size; 77 rank = c->rank; 78 Mbs = c->Mbs; 79 80 PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag1)); 81 PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag2)); 82 83 PetscCall(PetscMalloc2(imax, (PetscInt ***)&idx, imax, &n)); 84 85 for (i = 0; i < imax; i++) { 86 PetscCall(ISGetIndices(is[i], &idx[i])); 87 PetscCall(ISGetLocalSize(is[i], &n[i])); 88 } 89 90 /* evaluate communication - mesg to who,length of mesg, and buffer space 91 required. Based on this, buffers are allocated, and data copied into them*/ 92 PetscCall(PetscCalloc4(size, &w1, size, &w2, size, &w3, size, &w4)); 93 for (i = 0; i < imax; i++) { 94 PetscCall(PetscArrayzero(w4, size)); /* initialise work vector*/ 95 idx_i = idx[i]; 96 len = n[i]; 97 for (j = 0; j < len; j++) { 98 row = idx_i[j]; 99 PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index set cannot have negative entries"); 100 PetscCall(PetscLayoutFindOwner(C->rmap, row * C->rmap->bs, &proc)); 101 w4[proc]++; 102 } 103 for (j = 0; j < size; j++) { 104 if (w4[j]) { 105 w1[j] += w4[j]; 106 w3[j]++; 107 } 108 } 109 } 110 111 nrqs = 0; /* no of outgoing messages */ 112 msz = 0; /* total mesg length (for all proc */ 113 w1[rank] = 0; /* no mesg sent to itself */ 114 w3[rank] = 0; 115 for (i = 0; i < size; i++) { 116 if (w1[i]) { 117 w2[i] = 1; 118 nrqs++; 119 } /* there exists a message to proc i */ 120 } 121 /* pa - is list of processors to communicate with */ 122 PetscCall(PetscMalloc1(nrqs, &pa)); 123 for (i = 0, j = 0; i < size; i++) { 124 if (w1[i]) { 125 pa[j] = i; 126 j++; 127 } 128 } 129 130 /* Each message would have a header = 1 + 2*(no of IS) + data */ 131 for (i = 0; i < nrqs; i++) { 132 j = pa[i]; 133 w1[j] += w2[j] + 2 * w3[j]; 134 msz += w1[j]; 135 } 136 137 /* Determine the number of messages to expect, their lengths, from from-ids */ 138 PetscCall(PetscGatherNumberOfMessages(comm, w2, w1, &nrqr)); 139 PetscCall(PetscGatherMessageLengths(comm, nrqs, nrqr, w1, &onodes1, &olengths1)); 140 141 /* Now post the Irecvs corresponding to these messages */ 142 PetscCall(PetscPostIrecvInt(comm, tag1, nrqr, onodes1, olengths1, &rbuf, &r_waits1)); 143 144 /* Allocate Memory for outgoing messages */ 145 PetscCall(PetscMalloc4(size, &outdat, size, &ptr, msz, &tmp, size, &ctr)); 146 PetscCall(PetscArrayzero(outdat, size)); 147 PetscCall(PetscArrayzero(ptr, size)); 148 { 149 PetscInt *iptr = tmp, ict = 0; 150 for (i = 0; i < nrqs; i++) { 151 j = pa[i]; 152 iptr += ict; 153 outdat[j] = iptr; 154 ict = w1[j]; 155 } 156 } 157 158 /* Form the outgoing messages */ 159 /*plug in the headers*/ 160 for (i = 0; i < nrqs; i++) { 161 j = pa[i]; 162 outdat[j][0] = 0; 163 PetscCall(PetscArrayzero(outdat[j] + 1, 2 * w3[j])); 164 ptr[j] = outdat[j] + 2 * w3[j] + 1; 165 } 166 167 /* Memory for doing local proc's work*/ 168 { 169 PetscCall(PetscCalloc5(imax, &table, imax, &data, imax, &isz, Mbs * imax, &d_p, (Mbs / PETSC_BITS_PER_BYTE + 1) * imax, &t_p)); 170 171 for (i = 0; i < imax; i++) { 172 table[i] = t_p + (Mbs / PETSC_BITS_PER_BYTE + 1) * i; 173 data[i] = d_p + (Mbs)*i; 174 } 175 } 176 177 /* Parse the IS and update local tables and the outgoing buf with the data*/ 178 { 179 PetscInt n_i, *data_i, isz_i, *outdat_j, ctr_j; 180 PetscBT table_i; 181 182 for (i = 0; i < imax; i++) { 183 PetscCall(PetscArrayzero(ctr, size)); 184 n_i = n[i]; 185 table_i = table[i]; 186 idx_i = idx[i]; 187 data_i = data[i]; 188 isz_i = isz[i]; 189 for (j = 0; j < n_i; j++) { /* parse the indices of each IS */ 190 row = idx_i[j]; 191 PetscCall(PetscLayoutFindOwner(C->rmap, row * C->rmap->bs, &proc)); 192 if (proc != rank) { /* copy to the outgoing buffer */ 193 ctr[proc]++; 194 *ptr[proc] = row; 195 ptr[proc]++; 196 } else { /* Update the local table */ 197 if (!PetscBTLookupSet(table_i, row)) data_i[isz_i++] = row; 198 } 199 } 200 /* Update the headers for the current IS */ 201 for (j = 0; j < size; j++) { /* Can Optimise this loop by using pa[] */ 202 if ((ctr_j = ctr[j])) { 203 outdat_j = outdat[j]; 204 k = ++outdat_j[0]; 205 outdat_j[2 * k] = ctr_j; 206 outdat_j[2 * k - 1] = i; 207 } 208 } 209 isz[i] = isz_i; 210 } 211 } 212 213 /* Now post the sends */ 214 PetscCall(PetscMalloc1(nrqs, &s_waits1)); 215 for (i = 0; i < nrqs; ++i) { 216 j = pa[i]; 217 PetscCallMPI(MPI_Isend(outdat[j], w1[j], MPIU_INT, j, tag1, comm, s_waits1 + i)); 218 } 219 220 /* No longer need the original indices*/ 221 for (i = 0; i < imax; ++i) PetscCall(ISRestoreIndices(is[i], idx + i)); 222 PetscCall(PetscFree2(*(PetscInt ***)&idx, n)); 223 224 PetscCall(PetscMalloc1(imax, &iscomms)); 225 for (i = 0; i < imax; ++i) { 226 PetscCall(PetscCommDuplicate(PetscObjectComm((PetscObject)is[i]), &iscomms[i], NULL)); 227 PetscCall(ISDestroy(&is[i])); 228 } 229 230 /* Do Local work*/ 231 PetscCall(MatIncreaseOverlap_MPIBAIJ_Local(C, imax, table, isz, data)); 232 233 /* Receive messages*/ 234 PetscCallMPI(MPI_Waitall(nrqr, r_waits1, MPI_STATUSES_IGNORE)); 235 PetscCallMPI(MPI_Waitall(nrqs, s_waits1, MPI_STATUSES_IGNORE)); 236 237 /* Phase 1 sends are complete - deallocate buffers */ 238 PetscCall(PetscFree4(outdat, ptr, tmp, ctr)); 239 PetscCall(PetscFree4(w1, w2, w3, w4)); 240 241 PetscCall(PetscMalloc1(nrqr, &xdata)); 242 PetscCall(PetscMalloc1(nrqr, &isz1)); 243 PetscCall(MatIncreaseOverlap_MPIBAIJ_Receive(C, nrqr, rbuf, xdata, isz1)); 244 if (rbuf) { 245 PetscCall(PetscFree(rbuf[0])); 246 PetscCall(PetscFree(rbuf)); 247 } 248 249 /* Send the data back*/ 250 /* Do a global reduction to know the buffer space req for incoming messages*/ 251 { 252 PetscMPIInt *rw1; 253 254 PetscCall(PetscCalloc1(size, &rw1)); 255 256 for (i = 0; i < nrqr; ++i) { 257 proc = onodes1[i]; 258 rw1[proc] = isz1[i]; 259 } 260 261 /* Determine the number of messages to expect, their lengths, from from-ids */ 262 PetscCall(PetscGatherMessageLengths(comm, nrqr, nrqs, rw1, &onodes2, &olengths2)); 263 PetscCall(PetscFree(rw1)); 264 } 265 /* Now post the Irecvs corresponding to these messages */ 266 PetscCall(PetscPostIrecvInt(comm, tag2, nrqs, onodes2, olengths2, &rbuf2, &r_waits2)); 267 268 /* Now post the sends */ 269 PetscCall(PetscMalloc1(nrqr, &s_waits2)); 270 for (i = 0; i < nrqr; ++i) { 271 j = onodes1[i]; 272 PetscCallMPI(MPI_Isend(xdata[i], isz1[i], MPIU_INT, j, tag2, comm, s_waits2 + i)); 273 } 274 275 PetscCall(PetscFree(onodes1)); 276 PetscCall(PetscFree(olengths1)); 277 278 /* receive work done on other processors*/ 279 { 280 PetscMPIInt idex; 281 PetscInt is_no, ct1, max, *rbuf2_i, isz_i, *data_i, jmax; 282 PetscBT table_i; 283 284 for (i = 0; i < nrqs; ++i) { 285 PetscCallMPI(MPI_Waitany(nrqs, r_waits2, &idex, MPI_STATUS_IGNORE)); 286 /* Process the message*/ 287 rbuf2_i = rbuf2[idex]; 288 ct1 = 2 * rbuf2_i[0] + 1; 289 jmax = rbuf2[idex][0]; 290 for (j = 1; j <= jmax; j++) { 291 max = rbuf2_i[2 * j]; 292 is_no = rbuf2_i[2 * j - 1]; 293 isz_i = isz[is_no]; 294 data_i = data[is_no]; 295 table_i = table[is_no]; 296 for (k = 0; k < max; k++, ct1++) { 297 row = rbuf2_i[ct1]; 298 if (!PetscBTLookupSet(table_i, row)) data_i[isz_i++] = row; 299 } 300 isz[is_no] = isz_i; 301 } 302 } 303 PetscCallMPI(MPI_Waitall(nrqr, s_waits2, MPI_STATUSES_IGNORE)); 304 } 305 306 for (i = 0; i < imax; ++i) { 307 PetscCall(ISCreateGeneral(iscomms[i], isz[i], data[i], PETSC_COPY_VALUES, is + i)); 308 PetscCall(PetscCommDestroy(&iscomms[i])); 309 } 310 311 PetscCall(PetscFree(iscomms)); 312 PetscCall(PetscFree(onodes2)); 313 PetscCall(PetscFree(olengths2)); 314 315 PetscCall(PetscFree(pa)); 316 if (rbuf2) { 317 PetscCall(PetscFree(rbuf2[0])); 318 PetscCall(PetscFree(rbuf2)); 319 } 320 PetscCall(PetscFree(s_waits1)); 321 PetscCall(PetscFree(r_waits1)); 322 PetscCall(PetscFree(s_waits2)); 323 PetscCall(PetscFree(r_waits2)); 324 PetscCall(PetscFree5(table, data, isz, d_p, t_p)); 325 if (xdata) { 326 PetscCall(PetscFree(xdata[0])); 327 PetscCall(PetscFree(xdata)); 328 } 329 PetscCall(PetscFree(isz1)); 330 PetscFunctionReturn(PETSC_SUCCESS); 331 } 332 333 /* 334 MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 335 the work on the local processor. 336 337 Inputs: 338 C - MAT_MPIBAIJ; 339 imax - total no of index sets processed at a time; 340 table - an array of char - size = Mbs bits. 341 342 Output: 343 isz - array containing the count of the solution elements corresponding 344 to each index set; 345 data - pointer to the solutions 346 */ 347 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C, PetscInt imax, PetscBT *table, PetscInt *isz, PetscInt **data) 348 { 349 Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 350 Mat A = c->A, B = c->B; 351 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data; 352 PetscInt start, end, val, max, rstart, cstart, *ai, *aj; 353 PetscInt *bi, *bj, *garray, i, j, k, row, *data_i, isz_i; 354 PetscBT table_i; 355 356 PetscFunctionBegin; 357 rstart = c->rstartbs; 358 cstart = c->cstartbs; 359 ai = a->i; 360 aj = a->j; 361 bi = b->i; 362 bj = b->j; 363 garray = c->garray; 364 365 for (i = 0; i < imax; i++) { 366 data_i = data[i]; 367 table_i = table[i]; 368 isz_i = isz[i]; 369 for (j = 0, max = isz[i]; j < max; j++) { 370 row = data_i[j] - rstart; 371 start = ai[row]; 372 end = ai[row + 1]; 373 for (k = start; k < end; k++) { /* Amat */ 374 val = aj[k] + cstart; 375 if (!PetscBTLookupSet(table_i, val)) data_i[isz_i++] = val; 376 } 377 start = bi[row]; 378 end = bi[row + 1]; 379 for (k = start; k < end; k++) { /* Bmat */ 380 val = garray[bj[k]]; 381 if (!PetscBTLookupSet(table_i, val)) data_i[isz_i++] = val; 382 } 383 } 384 isz[i] = isz_i; 385 } 386 PetscFunctionReturn(PETSC_SUCCESS); 387 } 388 /* 389 MatIncreaseOverlap_MPIBAIJ_Receive - Process the received messages, 390 and return the output 391 392 Input: 393 C - the matrix 394 nrqr - no of messages being processed. 395 rbuf - an array of pointers to the received requests 396 397 Output: 398 xdata - array of messages to be sent back 399 isz1 - size of each message 400 401 For better efficiency perhaps we should malloc separately each xdata[i], 402 then if a remalloc is required we need only copy the data for that one row 403 rather than all previous rows as it is now where a single large chunk of 404 memory is used. 405 406 */ 407 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C, PetscInt nrqr, PetscInt **rbuf, PetscInt **xdata, PetscInt *isz1) 408 { 409 Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 410 Mat A = c->A, B = c->B; 411 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data; 412 PetscInt rstart, cstart, *ai, *aj, *bi, *bj, *garray, i, j, k; 413 PetscInt row, total_sz, ct, ct1, ct2, ct3, mem_estimate, oct2, l, start, end; 414 PetscInt val, max1, max2, Mbs, no_malloc = 0, *tmp, new_estimate, ctr; 415 PetscInt *rbuf_i, kmax, rbuf_0; 416 PetscBT xtable; 417 418 PetscFunctionBegin; 419 Mbs = c->Mbs; 420 rstart = c->rstartbs; 421 cstart = c->cstartbs; 422 ai = a->i; 423 aj = a->j; 424 bi = b->i; 425 bj = b->j; 426 garray = c->garray; 427 428 for (i = 0, ct = 0, total_sz = 0; i < nrqr; ++i) { 429 rbuf_i = rbuf[i]; 430 rbuf_0 = rbuf_i[0]; 431 ct += rbuf_0; 432 for (j = 1; j <= rbuf_0; j++) total_sz += rbuf_i[2 * j]; 433 } 434 435 if (c->Mbs) max1 = ct * (a->nz + b->nz) / c->Mbs; 436 else max1 = 1; 437 mem_estimate = 3 * ((total_sz > max1 ? total_sz : max1) + 1); 438 if (nrqr) { 439 PetscCall(PetscMalloc1(mem_estimate, &xdata[0])); 440 ++no_malloc; 441 } 442 PetscCall(PetscBTCreate(Mbs, &xtable)); 443 PetscCall(PetscArrayzero(isz1, nrqr)); 444 445 ct3 = 0; 446 for (i = 0; i < nrqr; i++) { /* for easch mesg from proc i */ 447 rbuf_i = rbuf[i]; 448 rbuf_0 = rbuf_i[0]; 449 ct1 = 2 * rbuf_0 + 1; 450 ct2 = ct1; 451 ct3 += ct1; 452 for (j = 1; j <= rbuf_0; j++) { /* for each IS from proc i*/ 453 PetscCall(PetscBTMemzero(Mbs, xtable)); 454 oct2 = ct2; 455 kmax = rbuf_i[2 * j]; 456 for (k = 0; k < kmax; k++, ct1++) { 457 row = rbuf_i[ct1]; 458 if (!PetscBTLookupSet(xtable, row)) { 459 if (!(ct3 < mem_estimate)) { 460 new_estimate = (PetscInt)(1.5 * mem_estimate) + 1; 461 PetscCall(PetscMalloc1(new_estimate, &tmp)); 462 PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate)); 463 PetscCall(PetscFree(xdata[0])); 464 xdata[0] = tmp; 465 mem_estimate = new_estimate; 466 ++no_malloc; 467 for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1]; 468 } 469 xdata[i][ct2++] = row; 470 ct3++; 471 } 472 } 473 for (k = oct2, max2 = ct2; k < max2; k++) { 474 row = xdata[i][k] - rstart; 475 start = ai[row]; 476 end = ai[row + 1]; 477 for (l = start; l < end; l++) { 478 val = aj[l] + cstart; 479 if (!PetscBTLookupSet(xtable, val)) { 480 if (!(ct3 < mem_estimate)) { 481 new_estimate = (PetscInt)(1.5 * mem_estimate) + 1; 482 PetscCall(PetscMalloc1(new_estimate, &tmp)); 483 PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate)); 484 PetscCall(PetscFree(xdata[0])); 485 xdata[0] = tmp; 486 mem_estimate = new_estimate; 487 ++no_malloc; 488 for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1]; 489 } 490 xdata[i][ct2++] = val; 491 ct3++; 492 } 493 } 494 start = bi[row]; 495 end = bi[row + 1]; 496 for (l = start; l < end; l++) { 497 val = garray[bj[l]]; 498 if (!PetscBTLookupSet(xtable, val)) { 499 if (!(ct3 < mem_estimate)) { 500 new_estimate = (PetscInt)(1.5 * mem_estimate) + 1; 501 PetscCall(PetscMalloc1(new_estimate, &tmp)); 502 PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate)); 503 PetscCall(PetscFree(xdata[0])); 504 xdata[0] = tmp; 505 mem_estimate = new_estimate; 506 ++no_malloc; 507 for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1]; 508 } 509 xdata[i][ct2++] = val; 510 ct3++; 511 } 512 } 513 } 514 /* Update the header*/ 515 xdata[i][2 * j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 516 xdata[i][2 * j - 1] = rbuf_i[2 * j - 1]; 517 } 518 xdata[i][0] = rbuf_0; 519 if (i + 1 < nrqr) xdata[i + 1] = xdata[i] + ct2; 520 isz1[i] = ct2; /* size of each message */ 521 } 522 PetscCall(PetscBTDestroy(&xtable)); 523 PetscCall(PetscInfo(C, "Allocated %" PetscInt_FMT " bytes, required %" PetscInt_FMT ", no of mallocs = %" PetscInt_FMT "\n", mem_estimate, ct3, no_malloc)); 524 PetscFunctionReturn(PETSC_SUCCESS); 525 } 526 527 PetscErrorCode MatCreateSubMatrices_MPIBAIJ(Mat C, PetscInt ismax, const IS isrow[], const IS iscol[], MatReuse scall, Mat *submat[]) 528 { 529 IS *isrow_block, *iscol_block; 530 Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 531 PetscInt nmax, nstages_local, nstages, i, pos, max_no, N = C->cmap->N, bs = C->rmap->bs; 532 Mat_SeqBAIJ *subc; 533 Mat_SubSppt *smat; 534 535 PetscFunctionBegin; 536 /* The compression and expansion should be avoided. Doesn't point 537 out errors, might change the indices, hence buggey */ 538 PetscCall(PetscMalloc2(ismax, &isrow_block, ismax, &iscol_block)); 539 PetscCall(ISCompressIndicesGeneral(C->rmap->N, C->rmap->n, bs, ismax, isrow, isrow_block)); 540 if (isrow == iscol) { 541 for (i = 0; i < ismax; i++) { 542 iscol_block[i] = isrow_block[i]; 543 PetscCall(PetscObjectReference((PetscObject)iscol_block[i])); 544 } 545 } else PetscCall(ISCompressIndicesGeneral(N, C->cmap->n, bs, ismax, iscol, iscol_block)); 546 547 /* Determine the number of stages through which submatrices are done */ 548 if (!C->cmap->N) nmax = 20 * 1000000 / sizeof(PetscInt); 549 else nmax = 20 * 1000000 / (c->Nbs * sizeof(PetscInt)); 550 if (!nmax) nmax = 1; 551 552 if (scall == MAT_INITIAL_MATRIX) { 553 nstages_local = ismax / nmax + ((ismax % nmax) ? 1 : 0); /* local nstages */ 554 555 /* Make sure every processor loops through the nstages */ 556 PetscCall(MPIU_Allreduce(&nstages_local, &nstages, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)C))); 557 558 /* Allocate memory to hold all the submatrices and dummy submatrices */ 559 PetscCall(PetscCalloc1(ismax + nstages, submat)); 560 } else { /* MAT_REUSE_MATRIX */ 561 if (ismax) { 562 subc = (Mat_SeqBAIJ *)((*submat)[0]->data); 563 smat = subc->submatis1; 564 } else { /* (*submat)[0] is a dummy matrix */ 565 smat = (Mat_SubSppt *)(*submat)[0]->data; 566 } 567 PetscCheck(smat, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "MatCreateSubMatrices(...,MAT_REUSE_MATRIX,...) requires submat"); 568 nstages = smat->nstages; 569 } 570 571 for (i = 0, pos = 0; i < nstages; i++) { 572 if (pos + nmax <= ismax) max_no = nmax; 573 else if (pos >= ismax) max_no = 0; 574 else max_no = ismax - pos; 575 576 PetscCall(MatCreateSubMatrices_MPIBAIJ_local(C, max_no, isrow_block + pos, iscol_block + pos, scall, *submat + pos)); 577 if (!max_no) { 578 if (scall == MAT_INITIAL_MATRIX) { /* submat[pos] is a dummy matrix */ 579 smat = (Mat_SubSppt *)(*submat)[pos]->data; 580 smat->nstages = nstages; 581 } 582 pos++; /* advance to next dummy matrix if any */ 583 } else pos += max_no; 584 } 585 586 if (scall == MAT_INITIAL_MATRIX && ismax) { 587 /* save nstages for reuse */ 588 subc = (Mat_SeqBAIJ *)((*submat)[0]->data); 589 smat = subc->submatis1; 590 smat->nstages = nstages; 591 } 592 593 for (i = 0; i < ismax; i++) { 594 PetscCall(ISDestroy(&isrow_block[i])); 595 PetscCall(ISDestroy(&iscol_block[i])); 596 } 597 PetscCall(PetscFree2(isrow_block, iscol_block)); 598 PetscFunctionReturn(PETSC_SUCCESS); 599 } 600 601 /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 602 PetscErrorCode MatCreateSubMatrices_MPIBAIJ_local(Mat C, PetscInt ismax, const IS isrow[], const IS iscol[], MatReuse scall, Mat *submats) 603 { 604 Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 605 Mat A = c->A; 606 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)c->B->data, *subc; 607 const PetscInt **icol, **irow; 608 PetscInt *nrow, *ncol, start; 609 PetscMPIInt rank, size, tag0, tag2, tag3, tag4, *w1, *w2, *w3, *w4, nrqr; 610 PetscInt **sbuf1, **sbuf2, *sbuf2_i, i, j, k, l, ct1, ct2, **rbuf1, row, proc = -1; 611 PetscInt nrqs = 0, msz, **ptr = NULL, *req_size = NULL, *ctr = NULL, *pa, *tmp = NULL, tcol; 612 PetscInt **rbuf3 = NULL, *req_source1 = NULL, *req_source2, **sbuf_aj, **rbuf2 = NULL, max1, max2; 613 PetscInt **lens, is_no, ncols, *cols, mat_i, *mat_j, tmp2, jmax; 614 #if defined(PETSC_USE_CTABLE) 615 PetscHMapI *cmap, cmap_i = NULL, *rmap, rmap_i; 616 #else 617 PetscInt **cmap, *cmap_i = NULL, **rmap, *rmap_i; 618 #endif 619 const PetscInt *irow_i, *icol_i; 620 PetscInt ctr_j, *sbuf1_j, *sbuf_aj_i, *rbuf1_i, kmax, *lens_i; 621 MPI_Request *s_waits1, *r_waits1, *s_waits2, *r_waits2, *r_waits3; 622 MPI_Request *r_waits4, *s_waits3, *s_waits4; 623 MPI_Comm comm; 624 PetscScalar **rbuf4, *rbuf4_i = NULL, **sbuf_aa, *vals, *mat_a = NULL, *imat_a = NULL, *sbuf_aa_i; 625 PetscMPIInt *onodes1, *olengths1, end; 626 PetscInt **row2proc, *row2proc_i, *imat_ilen, *imat_j, *imat_i; 627 Mat_SubSppt *smat_i; 628 PetscBool *issorted, colflag, iscsorted = PETSC_TRUE; 629 PetscInt *sbuf1_i, *rbuf2_i, *rbuf3_i, ilen; 630 PetscInt bs = C->rmap->bs, bs2 = c->bs2, rstart = c->rstartbs; 631 PetscBool ijonly = c->ijonly; /* private flag indicates only matrix data structures are requested */ 632 PetscInt nzA, nzB, *a_i = a->i, *b_i = b->i, *a_j = a->j, *b_j = b->j, ctmp, imark, *cworkA, *cworkB; 633 PetscScalar *vworkA = NULL, *vworkB = NULL, *a_a = a->a, *b_a = b->a; 634 PetscInt cstart = c->cstartbs, *bmap = c->garray; 635 PetscBool *allrows, *allcolumns; 636 637 PetscFunctionBegin; 638 PetscCall(PetscObjectGetComm((PetscObject)C, &comm)); 639 size = c->size; 640 rank = c->rank; 641 642 PetscCall(PetscMalloc5(ismax, &row2proc, ismax, &cmap, ismax, &rmap, ismax + 1, &allcolumns, ismax, &allrows)); 643 PetscCall(PetscMalloc5(ismax, (PetscInt ***)&irow, ismax, (PetscInt ***)&icol, ismax, &nrow, ismax, &ncol, ismax, &issorted)); 644 645 for (i = 0; i < ismax; i++) { 646 PetscCall(ISSorted(iscol[i], &issorted[i])); 647 if (!issorted[i]) iscsorted = issorted[i]; /* columns are not sorted! */ 648 PetscCall(ISSorted(isrow[i], &issorted[i])); 649 650 /* Check for special case: allcolumns */ 651 PetscCall(ISIdentity(iscol[i], &colflag)); 652 PetscCall(ISGetLocalSize(iscol[i], &ncol[i])); 653 654 if (colflag && ncol[i] == c->Nbs) { 655 allcolumns[i] = PETSC_TRUE; 656 icol[i] = NULL; 657 } else { 658 allcolumns[i] = PETSC_FALSE; 659 PetscCall(ISGetIndices(iscol[i], &icol[i])); 660 } 661 662 /* Check for special case: allrows */ 663 PetscCall(ISIdentity(isrow[i], &colflag)); 664 PetscCall(ISGetLocalSize(isrow[i], &nrow[i])); 665 if (colflag && nrow[i] == c->Mbs) { 666 allrows[i] = PETSC_TRUE; 667 irow[i] = NULL; 668 } else { 669 allrows[i] = PETSC_FALSE; 670 PetscCall(ISGetIndices(isrow[i], &irow[i])); 671 } 672 } 673 674 if (scall == MAT_REUSE_MATRIX) { 675 /* Assumes new rows are same length as the old rows */ 676 for (i = 0; i < ismax; i++) { 677 subc = (Mat_SeqBAIJ *)(submats[i]->data); 678 PetscCheck(subc->mbs == nrow[i] && subc->nbs == ncol[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong size"); 679 680 /* Initial matrix as if empty */ 681 PetscCall(PetscArrayzero(subc->ilen, subc->mbs)); 682 683 /* Initial matrix as if empty */ 684 submats[i]->factortype = C->factortype; 685 686 smat_i = subc->submatis1; 687 688 nrqs = smat_i->nrqs; 689 nrqr = smat_i->nrqr; 690 rbuf1 = smat_i->rbuf1; 691 rbuf2 = smat_i->rbuf2; 692 rbuf3 = smat_i->rbuf3; 693 req_source2 = smat_i->req_source2; 694 695 sbuf1 = smat_i->sbuf1; 696 sbuf2 = smat_i->sbuf2; 697 ptr = smat_i->ptr; 698 tmp = smat_i->tmp; 699 ctr = smat_i->ctr; 700 701 pa = smat_i->pa; 702 req_size = smat_i->req_size; 703 req_source1 = smat_i->req_source1; 704 705 allcolumns[i] = smat_i->allcolumns; 706 allrows[i] = smat_i->allrows; 707 row2proc[i] = smat_i->row2proc; 708 rmap[i] = smat_i->rmap; 709 cmap[i] = smat_i->cmap; 710 } 711 712 if (!ismax) { /* Get dummy submatrices and retrieve struct submatis1 */ 713 PetscCheck(submats[0], PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "submats are null, cannot reuse"); 714 smat_i = (Mat_SubSppt *)submats[0]->data; 715 716 nrqs = smat_i->nrqs; 717 nrqr = smat_i->nrqr; 718 rbuf1 = smat_i->rbuf1; 719 rbuf2 = smat_i->rbuf2; 720 rbuf3 = smat_i->rbuf3; 721 req_source2 = smat_i->req_source2; 722 723 sbuf1 = smat_i->sbuf1; 724 sbuf2 = smat_i->sbuf2; 725 ptr = smat_i->ptr; 726 tmp = smat_i->tmp; 727 ctr = smat_i->ctr; 728 729 pa = smat_i->pa; 730 req_size = smat_i->req_size; 731 req_source1 = smat_i->req_source1; 732 733 allcolumns[0] = PETSC_FALSE; 734 } 735 } else { /* scall == MAT_INITIAL_MATRIX */ 736 /* Get some new tags to keep the communication clean */ 737 PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag2)); 738 PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag3)); 739 740 /* evaluate communication - mesg to who, length of mesg, and buffer space 741 required. Based on this, buffers are allocated, and data copied into them*/ 742 PetscCall(PetscCalloc4(size, &w1, size, &w2, size, &w3, size, &w4)); /* mesg size, initialize work vectors */ 743 744 for (i = 0; i < ismax; i++) { 745 jmax = nrow[i]; 746 irow_i = irow[i]; 747 748 PetscCall(PetscMalloc1(jmax, &row2proc_i)); 749 row2proc[i] = row2proc_i; 750 751 if (issorted[i]) proc = 0; 752 for (j = 0; j < jmax; j++) { 753 if (!issorted[i]) proc = 0; 754 if (allrows[i]) row = j; 755 else row = irow_i[j]; 756 757 while (row >= c->rangebs[proc + 1]) proc++; 758 w4[proc]++; 759 row2proc_i[j] = proc; /* map row index to proc */ 760 } 761 for (j = 0; j < size; j++) { 762 if (w4[j]) { 763 w1[j] += w4[j]; 764 w3[j]++; 765 w4[j] = 0; 766 } 767 } 768 } 769 770 nrqs = 0; /* no of outgoing messages */ 771 msz = 0; /* total mesg length (for all procs) */ 772 w1[rank] = 0; /* no mesg sent to self */ 773 w3[rank] = 0; 774 for (i = 0; i < size; i++) { 775 if (w1[i]) { 776 w2[i] = 1; 777 nrqs++; 778 } /* there exists a message to proc i */ 779 } 780 PetscCall(PetscMalloc1(nrqs, &pa)); /*(proc -array)*/ 781 for (i = 0, j = 0; i < size; i++) { 782 if (w1[i]) { 783 pa[j] = i; 784 j++; 785 } 786 } 787 788 /* Each message would have a header = 1 + 2*(no of IS) + data */ 789 for (i = 0; i < nrqs; i++) { 790 j = pa[i]; 791 w1[j] += w2[j] + 2 * w3[j]; 792 msz += w1[j]; 793 } 794 PetscCall(PetscInfo(0, "Number of outgoing messages %" PetscInt_FMT " Total message length %" PetscInt_FMT "\n", nrqs, msz)); 795 796 /* Determine the number of messages to expect, their lengths, from from-ids */ 797 PetscCall(PetscGatherNumberOfMessages(comm, w2, w1, &nrqr)); 798 PetscCall(PetscGatherMessageLengths(comm, nrqs, nrqr, w1, &onodes1, &olengths1)); 799 800 /* Now post the Irecvs corresponding to these messages */ 801 PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag0)); 802 PetscCall(PetscPostIrecvInt(comm, tag0, nrqr, onodes1, olengths1, &rbuf1, &r_waits1)); 803 804 /* Allocate Memory for outgoing messages */ 805 PetscCall(PetscMalloc4(size, &sbuf1, size, &ptr, 2 * msz, &tmp, size, &ctr)); 806 PetscCall(PetscArrayzero(sbuf1, size)); 807 PetscCall(PetscArrayzero(ptr, size)); 808 809 { 810 PetscInt *iptr = tmp; 811 k = 0; 812 for (i = 0; i < nrqs; i++) { 813 j = pa[i]; 814 iptr += k; 815 sbuf1[j] = iptr; 816 k = w1[j]; 817 } 818 } 819 820 /* Form the outgoing messages. Initialize the header space */ 821 for (i = 0; i < nrqs; i++) { 822 j = pa[i]; 823 sbuf1[j][0] = 0; 824 PetscCall(PetscArrayzero(sbuf1[j] + 1, 2 * w3[j])); 825 ptr[j] = sbuf1[j] + 2 * w3[j] + 1; 826 } 827 828 /* Parse the isrow and copy data into outbuf */ 829 for (i = 0; i < ismax; i++) { 830 row2proc_i = row2proc[i]; 831 PetscCall(PetscArrayzero(ctr, size)); 832 irow_i = irow[i]; 833 jmax = nrow[i]; 834 for (j = 0; j < jmax; j++) { /* parse the indices of each IS */ 835 proc = row2proc_i[j]; 836 if (allrows[i]) row = j; 837 else row = irow_i[j]; 838 839 if (proc != rank) { /* copy to the outgoing buf*/ 840 ctr[proc]++; 841 *ptr[proc] = row; 842 ptr[proc]++; 843 } 844 } 845 /* Update the headers for the current IS */ 846 for (j = 0; j < size; j++) { /* Can Optimise this loop too */ 847 if ((ctr_j = ctr[j])) { 848 sbuf1_j = sbuf1[j]; 849 k = ++sbuf1_j[0]; 850 sbuf1_j[2 * k] = ctr_j; 851 sbuf1_j[2 * k - 1] = i; 852 } 853 } 854 } 855 856 /* Now post the sends */ 857 PetscCall(PetscMalloc1(nrqs, &s_waits1)); 858 for (i = 0; i < nrqs; ++i) { 859 j = pa[i]; 860 PetscCallMPI(MPI_Isend(sbuf1[j], w1[j], MPIU_INT, j, tag0, comm, s_waits1 + i)); 861 } 862 863 /* Post Receives to capture the buffer size */ 864 PetscCall(PetscMalloc1(nrqs, &r_waits2)); 865 PetscCall(PetscMalloc3(nrqs, &req_source2, nrqs, &rbuf2, nrqs, &rbuf3)); 866 if (nrqs) rbuf2[0] = tmp + msz; 867 for (i = 1; i < nrqs; ++i) rbuf2[i] = rbuf2[i - 1] + w1[pa[i - 1]]; 868 for (i = 0; i < nrqs; ++i) { 869 j = pa[i]; 870 PetscCallMPI(MPI_Irecv(rbuf2[i], w1[j], MPIU_INT, j, tag2, comm, r_waits2 + i)); 871 } 872 873 /* Send to other procs the buf size they should allocate */ 874 /* Receive messages*/ 875 PetscCall(PetscMalloc1(nrqr, &s_waits2)); 876 PetscCall(PetscMalloc3(nrqr, &sbuf2, nrqr, &req_size, nrqr, &req_source1)); 877 878 PetscCallMPI(MPI_Waitall(nrqr, r_waits1, MPI_STATUSES_IGNORE)); 879 for (i = 0; i < nrqr; ++i) { 880 req_size[i] = 0; 881 rbuf1_i = rbuf1[i]; 882 start = 2 * rbuf1_i[0] + 1; 883 end = olengths1[i]; 884 PetscCall(PetscMalloc1(end, &sbuf2[i])); 885 sbuf2_i = sbuf2[i]; 886 for (j = start; j < end; j++) { 887 row = rbuf1_i[j] - rstart; 888 ncols = a_i[row + 1] - a_i[row] + b_i[row + 1] - b_i[row]; 889 sbuf2_i[j] = ncols; 890 req_size[i] += ncols; 891 } 892 req_source1[i] = onodes1[i]; 893 /* form the header */ 894 sbuf2_i[0] = req_size[i]; 895 for (j = 1; j < start; j++) sbuf2_i[j] = rbuf1_i[j]; 896 897 PetscCallMPI(MPI_Isend(sbuf2_i, end, MPIU_INT, req_source1[i], tag2, comm, s_waits2 + i)); 898 } 899 900 PetscCall(PetscFree(onodes1)); 901 PetscCall(PetscFree(olengths1)); 902 903 PetscCall(PetscFree(r_waits1)); 904 PetscCall(PetscFree4(w1, w2, w3, w4)); 905 906 /* Receive messages*/ 907 PetscCall(PetscMalloc1(nrqs, &r_waits3)); 908 909 PetscCallMPI(MPI_Waitall(nrqs, r_waits2, MPI_STATUSES_IGNORE)); 910 for (i = 0; i < nrqs; ++i) { 911 PetscCall(PetscMalloc1(rbuf2[i][0], &rbuf3[i])); 912 req_source2[i] = pa[i]; 913 PetscCallMPI(MPI_Irecv(rbuf3[i], rbuf2[i][0], MPIU_INT, req_source2[i], tag3, comm, r_waits3 + i)); 914 } 915 PetscCall(PetscFree(r_waits2)); 916 917 /* Wait on sends1 and sends2 */ 918 PetscCallMPI(MPI_Waitall(nrqs, s_waits1, MPI_STATUSES_IGNORE)); 919 PetscCallMPI(MPI_Waitall(nrqr, s_waits2, MPI_STATUSES_IGNORE)); 920 PetscCall(PetscFree(s_waits1)); 921 PetscCall(PetscFree(s_waits2)); 922 923 /* Now allocate sending buffers for a->j, and send them off */ 924 PetscCall(PetscMalloc1(nrqr, &sbuf_aj)); 925 for (i = 0, j = 0; i < nrqr; i++) j += req_size[i]; 926 if (nrqr) PetscCall(PetscMalloc1(j, &sbuf_aj[0])); 927 for (i = 1; i < nrqr; i++) sbuf_aj[i] = sbuf_aj[i - 1] + req_size[i - 1]; 928 929 PetscCall(PetscMalloc1(nrqr, &s_waits3)); 930 { 931 for (i = 0; i < nrqr; i++) { 932 rbuf1_i = rbuf1[i]; 933 sbuf_aj_i = sbuf_aj[i]; 934 ct1 = 2 * rbuf1_i[0] + 1; 935 ct2 = 0; 936 for (j = 1, max1 = rbuf1_i[0]; j <= max1; j++) { 937 kmax = rbuf1[i][2 * j]; 938 for (k = 0; k < kmax; k++, ct1++) { 939 row = rbuf1_i[ct1] - rstart; 940 nzA = a_i[row + 1] - a_i[row]; 941 nzB = b_i[row + 1] - b_i[row]; 942 ncols = nzA + nzB; 943 cworkA = a_j + a_i[row]; 944 cworkB = b_j + b_i[row]; 945 946 /* load the column indices for this row into cols */ 947 cols = sbuf_aj_i + ct2; 948 for (l = 0; l < nzB; l++) { 949 if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 950 else break; 951 } 952 imark = l; 953 for (l = 0; l < nzA; l++) cols[imark + l] = cstart + cworkA[l]; 954 for (l = imark; l < nzB; l++) cols[nzA + l] = bmap[cworkB[l]]; 955 ct2 += ncols; 956 } 957 } 958 PetscCallMPI(MPI_Isend(sbuf_aj_i, req_size[i], MPIU_INT, req_source1[i], tag3, comm, s_waits3 + i)); 959 } 960 } 961 962 /* create col map: global col of C -> local col of submatrices */ 963 #if defined(PETSC_USE_CTABLE) 964 for (i = 0; i < ismax; i++) { 965 if (!allcolumns[i]) { 966 PetscCall(PetscHMapICreateWithSize(ncol[i], cmap + i)); 967 968 jmax = ncol[i]; 969 icol_i = icol[i]; 970 cmap_i = cmap[i]; 971 for (j = 0; j < jmax; j++) PetscCall(PetscHMapISet(cmap[i], icol_i[j] + 1, j + 1)); 972 } else cmap[i] = NULL; 973 } 974 #else 975 for (i = 0; i < ismax; i++) { 976 if (!allcolumns[i]) { 977 PetscCall(PetscCalloc1(c->Nbs, &cmap[i])); 978 jmax = ncol[i]; 979 icol_i = icol[i]; 980 cmap_i = cmap[i]; 981 for (j = 0; j < jmax; j++) cmap_i[icol_i[j]] = j + 1; 982 } else cmap[i] = NULL; 983 } 984 #endif 985 986 /* Create lens which is required for MatCreate... */ 987 for (i = 0, j = 0; i < ismax; i++) j += nrow[i]; 988 PetscCall(PetscMalloc1(ismax, &lens)); 989 990 if (ismax) PetscCall(PetscCalloc1(j, &lens[0])); 991 for (i = 1; i < ismax; i++) lens[i] = lens[i - 1] + nrow[i - 1]; 992 993 /* Update lens from local data */ 994 for (i = 0; i < ismax; i++) { 995 row2proc_i = row2proc[i]; 996 jmax = nrow[i]; 997 if (!allcolumns[i]) cmap_i = cmap[i]; 998 irow_i = irow[i]; 999 lens_i = lens[i]; 1000 for (j = 0; j < jmax; j++) { 1001 if (allrows[i]) row = j; 1002 else row = irow_i[j]; /* global blocked row of C */ 1003 1004 proc = row2proc_i[j]; 1005 if (proc == rank) { 1006 /* Get indices from matA and then from matB */ 1007 #if defined(PETSC_USE_CTABLE) 1008 PetscInt tt; 1009 #endif 1010 row = row - rstart; 1011 nzA = a_i[row + 1] - a_i[row]; 1012 nzB = b_i[row + 1] - b_i[row]; 1013 cworkA = a_j + a_i[row]; 1014 cworkB = b_j + b_i[row]; 1015 1016 if (!allcolumns[i]) { 1017 #if defined(PETSC_USE_CTABLE) 1018 for (k = 0; k < nzA; k++) { 1019 PetscCall(PetscHMapIGetWithDefault(cmap_i, cstart + cworkA[k] + 1, 0, &tt)); 1020 if (tt) lens_i[j]++; 1021 } 1022 for (k = 0; k < nzB; k++) { 1023 PetscCall(PetscHMapIGetWithDefault(cmap_i, bmap[cworkB[k]] + 1, 0, &tt)); 1024 if (tt) lens_i[j]++; 1025 } 1026 1027 #else 1028 for (k = 0; k < nzA; k++) { 1029 if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1030 } 1031 for (k = 0; k < nzB; k++) { 1032 if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1033 } 1034 #endif 1035 } else { /* allcolumns */ 1036 lens_i[j] = nzA + nzB; 1037 } 1038 } 1039 } 1040 } 1041 1042 /* Create row map: global row of C -> local row of submatrices */ 1043 for (i = 0; i < ismax; i++) { 1044 if (!allrows[i]) { 1045 #if defined(PETSC_USE_CTABLE) 1046 PetscCall(PetscHMapICreateWithSize(nrow[i], rmap + i)); 1047 irow_i = irow[i]; 1048 jmax = nrow[i]; 1049 for (j = 0; j < jmax; j++) { 1050 if (allrows[i]) { 1051 PetscCall(PetscHMapISet(rmap[i], j + 1, j + 1)); 1052 } else { 1053 PetscCall(PetscHMapISet(rmap[i], irow_i[j] + 1, j + 1)); 1054 } 1055 } 1056 #else 1057 PetscCall(PetscCalloc1(c->Mbs, &rmap[i])); 1058 rmap_i = rmap[i]; 1059 irow_i = irow[i]; 1060 jmax = nrow[i]; 1061 for (j = 0; j < jmax; j++) { 1062 if (allrows[i]) rmap_i[j] = j; 1063 else rmap_i[irow_i[j]] = j; 1064 } 1065 #endif 1066 } else rmap[i] = NULL; 1067 } 1068 1069 /* Update lens from offproc data */ 1070 { 1071 PetscInt *rbuf2_i, *rbuf3_i, *sbuf1_i; 1072 1073 PetscCallMPI(MPI_Waitall(nrqs, r_waits3, MPI_STATUSES_IGNORE)); 1074 for (tmp2 = 0; tmp2 < nrqs; tmp2++) { 1075 sbuf1_i = sbuf1[pa[tmp2]]; 1076 jmax = sbuf1_i[0]; 1077 ct1 = 2 * jmax + 1; 1078 ct2 = 0; 1079 rbuf2_i = rbuf2[tmp2]; 1080 rbuf3_i = rbuf3[tmp2]; 1081 for (j = 1; j <= jmax; j++) { 1082 is_no = sbuf1_i[2 * j - 1]; 1083 max1 = sbuf1_i[2 * j]; 1084 lens_i = lens[is_no]; 1085 if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1086 rmap_i = rmap[is_no]; 1087 for (k = 0; k < max1; k++, ct1++) { 1088 if (allrows[is_no]) { 1089 row = sbuf1_i[ct1]; 1090 } else { 1091 #if defined(PETSC_USE_CTABLE) 1092 PetscCall(PetscHMapIGetWithDefault(rmap_i, sbuf1_i[ct1] + 1, 0, &row)); 1093 row--; 1094 PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table"); 1095 #else 1096 row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1097 #endif 1098 } 1099 max2 = rbuf2_i[ct1]; 1100 for (l = 0; l < max2; l++, ct2++) { 1101 if (!allcolumns[is_no]) { 1102 #if defined(PETSC_USE_CTABLE) 1103 PetscCall(PetscHMapIGetWithDefault(cmap_i, rbuf3_i[ct2] + 1, 0, &tcol)); 1104 #else 1105 tcol = cmap_i[rbuf3_i[ct2]]; 1106 #endif 1107 if (tcol) lens_i[row]++; 1108 } else { /* allcolumns */ 1109 lens_i[row]++; /* lens_i[row] += max2 ? */ 1110 } 1111 } 1112 } 1113 } 1114 } 1115 } 1116 PetscCall(PetscFree(r_waits3)); 1117 PetscCallMPI(MPI_Waitall(nrqr, s_waits3, MPI_STATUSES_IGNORE)); 1118 PetscCall(PetscFree(s_waits3)); 1119 1120 /* Create the submatrices */ 1121 for (i = 0; i < ismax; i++) { 1122 PetscInt bs_tmp; 1123 if (ijonly) bs_tmp = 1; 1124 else bs_tmp = bs; 1125 1126 PetscCall(MatCreate(PETSC_COMM_SELF, submats + i)); 1127 PetscCall(MatSetSizes(submats[i], nrow[i] * bs_tmp, ncol[i] * bs_tmp, PETSC_DETERMINE, PETSC_DETERMINE)); 1128 1129 PetscCall(MatSetType(submats[i], ((PetscObject)A)->type_name)); 1130 PetscCall(MatSeqBAIJSetPreallocation(submats[i], bs_tmp, 0, lens[i])); 1131 PetscCall(MatSeqSBAIJSetPreallocation(submats[i], bs_tmp, 0, lens[i])); /* this subroutine is used by SBAIJ routines */ 1132 1133 /* create struct Mat_SubSppt and attached it to submat */ 1134 PetscCall(PetscNew(&smat_i)); 1135 subc = (Mat_SeqBAIJ *)submats[i]->data; 1136 subc->submatis1 = smat_i; 1137 1138 smat_i->destroy = submats[i]->ops->destroy; 1139 submats[i]->ops->destroy = MatDestroySubMatrix_SeqBAIJ; 1140 submats[i]->factortype = C->factortype; 1141 1142 smat_i->id = i; 1143 smat_i->nrqs = nrqs; 1144 smat_i->nrqr = nrqr; 1145 smat_i->rbuf1 = rbuf1; 1146 smat_i->rbuf2 = rbuf2; 1147 smat_i->rbuf3 = rbuf3; 1148 smat_i->sbuf2 = sbuf2; 1149 smat_i->req_source2 = req_source2; 1150 1151 smat_i->sbuf1 = sbuf1; 1152 smat_i->ptr = ptr; 1153 smat_i->tmp = tmp; 1154 smat_i->ctr = ctr; 1155 1156 smat_i->pa = pa; 1157 smat_i->req_size = req_size; 1158 smat_i->req_source1 = req_source1; 1159 1160 smat_i->allcolumns = allcolumns[i]; 1161 smat_i->allrows = allrows[i]; 1162 smat_i->singleis = PETSC_FALSE; 1163 smat_i->row2proc = row2proc[i]; 1164 smat_i->rmap = rmap[i]; 1165 smat_i->cmap = cmap[i]; 1166 } 1167 1168 if (!ismax) { /* Create dummy submats[0] for reuse struct subc */ 1169 PetscCall(MatCreate(PETSC_COMM_SELF, &submats[0])); 1170 PetscCall(MatSetSizes(submats[0], 0, 0, PETSC_DETERMINE, PETSC_DETERMINE)); 1171 PetscCall(MatSetType(submats[0], MATDUMMY)); 1172 1173 /* create struct Mat_SubSppt and attached it to submat */ 1174 PetscCall(PetscNew(&smat_i)); 1175 submats[0]->data = (void *)smat_i; 1176 1177 smat_i->destroy = submats[0]->ops->destroy; 1178 submats[0]->ops->destroy = MatDestroySubMatrix_Dummy; 1179 submats[0]->factortype = C->factortype; 1180 1181 smat_i->id = 0; 1182 smat_i->nrqs = nrqs; 1183 smat_i->nrqr = nrqr; 1184 smat_i->rbuf1 = rbuf1; 1185 smat_i->rbuf2 = rbuf2; 1186 smat_i->rbuf3 = rbuf3; 1187 smat_i->sbuf2 = sbuf2; 1188 smat_i->req_source2 = req_source2; 1189 1190 smat_i->sbuf1 = sbuf1; 1191 smat_i->ptr = ptr; 1192 smat_i->tmp = tmp; 1193 smat_i->ctr = ctr; 1194 1195 smat_i->pa = pa; 1196 smat_i->req_size = req_size; 1197 smat_i->req_source1 = req_source1; 1198 1199 smat_i->allcolumns = PETSC_FALSE; 1200 smat_i->singleis = PETSC_FALSE; 1201 smat_i->row2proc = NULL; 1202 smat_i->rmap = NULL; 1203 smat_i->cmap = NULL; 1204 } 1205 1206 if (ismax) PetscCall(PetscFree(lens[0])); 1207 PetscCall(PetscFree(lens)); 1208 if (sbuf_aj) { 1209 PetscCall(PetscFree(sbuf_aj[0])); 1210 PetscCall(PetscFree(sbuf_aj)); 1211 } 1212 1213 } /* endof scall == MAT_INITIAL_MATRIX */ 1214 1215 /* Post recv matrix values */ 1216 if (!ijonly) { 1217 PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag4)); 1218 PetscCall(PetscMalloc1(nrqs, &rbuf4)); 1219 PetscCall(PetscMalloc1(nrqs, &r_waits4)); 1220 for (i = 0; i < nrqs; ++i) { 1221 PetscCall(PetscMalloc1(rbuf2[i][0] * bs2, &rbuf4[i])); 1222 PetscCallMPI(MPI_Irecv(rbuf4[i], rbuf2[i][0] * bs2, MPIU_SCALAR, req_source2[i], tag4, comm, r_waits4 + i)); 1223 } 1224 1225 /* Allocate sending buffers for a->a, and send them off */ 1226 PetscCall(PetscMalloc1(nrqr, &sbuf_aa)); 1227 for (i = 0, j = 0; i < nrqr; i++) j += req_size[i]; 1228 1229 if (nrqr) PetscCall(PetscMalloc1(j * bs2, &sbuf_aa[0])); 1230 for (i = 1; i < nrqr; i++) sbuf_aa[i] = sbuf_aa[i - 1] + req_size[i - 1] * bs2; 1231 1232 PetscCall(PetscMalloc1(nrqr, &s_waits4)); 1233 1234 for (i = 0; i < nrqr; i++) { 1235 rbuf1_i = rbuf1[i]; 1236 sbuf_aa_i = sbuf_aa[i]; 1237 ct1 = 2 * rbuf1_i[0] + 1; 1238 ct2 = 0; 1239 for (j = 1, max1 = rbuf1_i[0]; j <= max1; j++) { 1240 kmax = rbuf1_i[2 * j]; 1241 for (k = 0; k < kmax; k++, ct1++) { 1242 row = rbuf1_i[ct1] - rstart; 1243 nzA = a_i[row + 1] - a_i[row]; 1244 nzB = b_i[row + 1] - b_i[row]; 1245 ncols = nzA + nzB; 1246 cworkB = b_j + b_i[row]; 1247 vworkA = a_a + a_i[row] * bs2; 1248 vworkB = b_a + b_i[row] * bs2; 1249 1250 /* load the column values for this row into vals*/ 1251 vals = sbuf_aa_i + ct2 * bs2; 1252 for (l = 0; l < nzB; l++) { 1253 if ((bmap[cworkB[l]]) < cstart) { 1254 PetscCall(PetscArraycpy(vals + l * bs2, vworkB + l * bs2, bs2)); 1255 } else break; 1256 } 1257 imark = l; 1258 for (l = 0; l < nzA; l++) PetscCall(PetscArraycpy(vals + (imark + l) * bs2, vworkA + l * bs2, bs2)); 1259 for (l = imark; l < nzB; l++) PetscCall(PetscArraycpy(vals + (nzA + l) * bs2, vworkB + l * bs2, bs2)); 1260 1261 ct2 += ncols; 1262 } 1263 } 1264 PetscCallMPI(MPI_Isend(sbuf_aa_i, req_size[i] * bs2, MPIU_SCALAR, req_source1[i], tag4, comm, s_waits4 + i)); 1265 } 1266 } 1267 1268 /* Assemble the matrices */ 1269 /* First assemble the local rows */ 1270 for (i = 0; i < ismax; i++) { 1271 row2proc_i = row2proc[i]; 1272 subc = (Mat_SeqBAIJ *)submats[i]->data; 1273 imat_ilen = subc->ilen; 1274 imat_j = subc->j; 1275 imat_i = subc->i; 1276 imat_a = subc->a; 1277 1278 if (!allcolumns[i]) cmap_i = cmap[i]; 1279 rmap_i = rmap[i]; 1280 irow_i = irow[i]; 1281 jmax = nrow[i]; 1282 for (j = 0; j < jmax; j++) { 1283 if (allrows[i]) row = j; 1284 else row = irow_i[j]; 1285 proc = row2proc_i[j]; 1286 1287 if (proc == rank) { 1288 row = row - rstart; 1289 nzA = a_i[row + 1] - a_i[row]; 1290 nzB = b_i[row + 1] - b_i[row]; 1291 cworkA = a_j + a_i[row]; 1292 cworkB = b_j + b_i[row]; 1293 if (!ijonly) { 1294 vworkA = a_a + a_i[row] * bs2; 1295 vworkB = b_a + b_i[row] * bs2; 1296 } 1297 1298 if (allrows[i]) { 1299 row = row + rstart; 1300 } else { 1301 #if defined(PETSC_USE_CTABLE) 1302 PetscCall(PetscHMapIGetWithDefault(rmap_i, row + rstart + 1, 0, &row)); 1303 row--; 1304 1305 PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table"); 1306 #else 1307 row = rmap_i[row + rstart]; 1308 #endif 1309 } 1310 mat_i = imat_i[row]; 1311 if (!ijonly) mat_a = imat_a + mat_i * bs2; 1312 mat_j = imat_j + mat_i; 1313 ilen = imat_ilen[row]; 1314 1315 /* load the column indices for this row into cols*/ 1316 if (!allcolumns[i]) { 1317 for (l = 0; l < nzB; l++) { 1318 if ((ctmp = bmap[cworkB[l]]) < cstart) { 1319 #if defined(PETSC_USE_CTABLE) 1320 PetscCall(PetscHMapIGetWithDefault(cmap_i, ctmp + 1, 0, &tcol)); 1321 if (tcol) { 1322 #else 1323 if ((tcol = cmap_i[ctmp])) { 1324 #endif 1325 *mat_j++ = tcol - 1; 1326 PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2)); 1327 mat_a += bs2; 1328 ilen++; 1329 } 1330 } else break; 1331 } 1332 imark = l; 1333 for (l = 0; l < nzA; l++) { 1334 #if defined(PETSC_USE_CTABLE) 1335 PetscCall(PetscHMapIGetWithDefault(cmap_i, cstart + cworkA[l] + 1, 0, &tcol)); 1336 if (tcol) { 1337 #else 1338 if ((tcol = cmap_i[cstart + cworkA[l]])) { 1339 #endif 1340 *mat_j++ = tcol - 1; 1341 if (!ijonly) { 1342 PetscCall(PetscArraycpy(mat_a, vworkA + l * bs2, bs2)); 1343 mat_a += bs2; 1344 } 1345 ilen++; 1346 } 1347 } 1348 for (l = imark; l < nzB; l++) { 1349 #if defined(PETSC_USE_CTABLE) 1350 PetscCall(PetscHMapIGetWithDefault(cmap_i, bmap[cworkB[l]] + 1, 0, &tcol)); 1351 if (tcol) { 1352 #else 1353 if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1354 #endif 1355 *mat_j++ = tcol - 1; 1356 if (!ijonly) { 1357 PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2)); 1358 mat_a += bs2; 1359 } 1360 ilen++; 1361 } 1362 } 1363 } else { /* allcolumns */ 1364 for (l = 0; l < nzB; l++) { 1365 if ((ctmp = bmap[cworkB[l]]) < cstart) { 1366 *mat_j++ = ctmp; 1367 PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2)); 1368 mat_a += bs2; 1369 ilen++; 1370 } else break; 1371 } 1372 imark = l; 1373 for (l = 0; l < nzA; l++) { 1374 *mat_j++ = cstart + cworkA[l]; 1375 if (!ijonly) { 1376 PetscCall(PetscArraycpy(mat_a, vworkA + l * bs2, bs2)); 1377 mat_a += bs2; 1378 } 1379 ilen++; 1380 } 1381 for (l = imark; l < nzB; l++) { 1382 *mat_j++ = bmap[cworkB[l]]; 1383 if (!ijonly) { 1384 PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2)); 1385 mat_a += bs2; 1386 } 1387 ilen++; 1388 } 1389 } 1390 imat_ilen[row] = ilen; 1391 } 1392 } 1393 } 1394 1395 /* Now assemble the off proc rows */ 1396 if (!ijonly) PetscCallMPI(MPI_Waitall(nrqs, r_waits4, MPI_STATUSES_IGNORE)); 1397 for (tmp2 = 0; tmp2 < nrqs; tmp2++) { 1398 sbuf1_i = sbuf1[pa[tmp2]]; 1399 jmax = sbuf1_i[0]; 1400 ct1 = 2 * jmax + 1; 1401 ct2 = 0; 1402 rbuf2_i = rbuf2[tmp2]; 1403 rbuf3_i = rbuf3[tmp2]; 1404 if (!ijonly) rbuf4_i = rbuf4[tmp2]; 1405 for (j = 1; j <= jmax; j++) { 1406 is_no = sbuf1_i[2 * j - 1]; 1407 rmap_i = rmap[is_no]; 1408 if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1409 subc = (Mat_SeqBAIJ *)submats[is_no]->data; 1410 imat_ilen = subc->ilen; 1411 imat_j = subc->j; 1412 imat_i = subc->i; 1413 if (!ijonly) imat_a = subc->a; 1414 max1 = sbuf1_i[2 * j]; 1415 for (k = 0; k < max1; k++, ct1++) { /* for each recved block row */ 1416 row = sbuf1_i[ct1]; 1417 1418 if (allrows[is_no]) { 1419 row = sbuf1_i[ct1]; 1420 } else { 1421 #if defined(PETSC_USE_CTABLE) 1422 PetscCall(PetscHMapIGetWithDefault(rmap_i, row + 1, 0, &row)); 1423 row--; 1424 PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table"); 1425 #else 1426 row = rmap_i[row]; 1427 #endif 1428 } 1429 ilen = imat_ilen[row]; 1430 mat_i = imat_i[row]; 1431 if (!ijonly) mat_a = imat_a + mat_i * bs2; 1432 mat_j = imat_j + mat_i; 1433 max2 = rbuf2_i[ct1]; 1434 if (!allcolumns[is_no]) { 1435 for (l = 0; l < max2; l++, ct2++) { 1436 #if defined(PETSC_USE_CTABLE) 1437 PetscCall(PetscHMapIGetWithDefault(cmap_i, rbuf3_i[ct2] + 1, 0, &tcol)); 1438 #else 1439 tcol = cmap_i[rbuf3_i[ct2]]; 1440 #endif 1441 if (tcol) { 1442 *mat_j++ = tcol - 1; 1443 if (!ijonly) { 1444 PetscCall(PetscArraycpy(mat_a, rbuf4_i + ct2 * bs2, bs2)); 1445 mat_a += bs2; 1446 } 1447 ilen++; 1448 } 1449 } 1450 } else { /* allcolumns */ 1451 for (l = 0; l < max2; l++, ct2++) { 1452 *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */ 1453 if (!ijonly) { 1454 PetscCall(PetscArraycpy(mat_a, rbuf4_i + ct2 * bs2, bs2)); 1455 mat_a += bs2; 1456 } 1457 ilen++; 1458 } 1459 } 1460 imat_ilen[row] = ilen; 1461 } 1462 } 1463 } 1464 1465 if (!iscsorted) { /* sort column indices of the rows */ 1466 MatScalar *work; 1467 1468 PetscCall(PetscMalloc1(bs2, &work)); 1469 for (i = 0; i < ismax; i++) { 1470 subc = (Mat_SeqBAIJ *)submats[i]->data; 1471 imat_ilen = subc->ilen; 1472 imat_j = subc->j; 1473 imat_i = subc->i; 1474 if (!ijonly) imat_a = subc->a; 1475 if (allcolumns[i]) continue; 1476 1477 jmax = nrow[i]; 1478 for (j = 0; j < jmax; j++) { 1479 mat_i = imat_i[j]; 1480 mat_j = imat_j + mat_i; 1481 ilen = imat_ilen[j]; 1482 if (ijonly) { 1483 PetscCall(PetscSortInt(ilen, mat_j)); 1484 } else { 1485 mat_a = imat_a + mat_i * bs2; 1486 PetscCall(PetscSortIntWithDataArray(ilen, mat_j, mat_a, bs2 * sizeof(MatScalar), work)); 1487 } 1488 } 1489 } 1490 PetscCall(PetscFree(work)); 1491 } 1492 1493 if (!ijonly) { 1494 PetscCall(PetscFree(r_waits4)); 1495 PetscCallMPI(MPI_Waitall(nrqr, s_waits4, MPI_STATUSES_IGNORE)); 1496 PetscCall(PetscFree(s_waits4)); 1497 } 1498 1499 /* Restore the indices */ 1500 for (i = 0; i < ismax; i++) { 1501 if (!allrows[i]) PetscCall(ISRestoreIndices(isrow[i], irow + i)); 1502 if (!allcolumns[i]) PetscCall(ISRestoreIndices(iscol[i], icol + i)); 1503 } 1504 1505 for (i = 0; i < ismax; i++) { 1506 PetscCall(MatAssemblyBegin(submats[i], MAT_FINAL_ASSEMBLY)); 1507 PetscCall(MatAssemblyEnd(submats[i], MAT_FINAL_ASSEMBLY)); 1508 } 1509 1510 PetscCall(PetscFree5(*(PetscInt ***)&irow, *(PetscInt ***)&icol, nrow, ncol, issorted)); 1511 PetscCall(PetscFree5(row2proc, cmap, rmap, allcolumns, allrows)); 1512 1513 if (!ijonly) { 1514 if (sbuf_aa) { 1515 PetscCall(PetscFree(sbuf_aa[0])); 1516 PetscCall(PetscFree(sbuf_aa)); 1517 } 1518 1519 for (i = 0; i < nrqs; ++i) PetscCall(PetscFree(rbuf4[i])); 1520 PetscCall(PetscFree(rbuf4)); 1521 } 1522 c->ijonly = PETSC_FALSE; /* set back to the default */ 1523 PetscFunctionReturn(PETSC_SUCCESS); 1524 } 1525