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 PetscErrorCode ierr; 17 PetscInt i,N=C->cmap->N, bs=C->rmap->bs; 18 IS *is_new; 19 20 PetscFunctionBegin; 21 ierr = PetscMalloc1(imax,&is_new);CHKERRQ(ierr); 22 /* Convert the indices into block format */ 23 ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,imax,is,is_new);CHKERRQ(ierr); 24 if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n"); 25 for (i=0; i<ov; ++i) { 26 ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr); 27 } 28 for (i=0; i<imax; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);} 29 ierr = ISExpandIndicesGeneral(N,N,bs,imax,is_new,is);CHKERRQ(ierr); 30 for (i=0; i<imax; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);} 31 ierr = PetscFree(is_new);CHKERRQ(ierr); 32 PetscFunctionReturn(0); 33 } 34 35 /* 36 Sample message format: 37 If a processor A wants processor B to process some elements corresponding 38 to index sets is[1], is[5] 39 mesg [0] = 2 (no of index sets in the mesg) 40 ----------- 41 mesg [1] = 1 => is[1] 42 mesg [2] = sizeof(is[1]); 43 ----------- 44 mesg [5] = 5 => is[5] 45 mesg [6] = sizeof(is[5]); 46 ----------- 47 mesg [7] 48 mesg [n] data(is[1]) 49 ----------- 50 mesg[n+1] 51 mesg[m] data(is[5]) 52 ----------- 53 54 Notes: 55 nrqs - no of requests sent (or to be sent out) 56 nrqr - no of requests recieved (which have to be or which have been processed 57 */ 58 PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[]) 59 { 60 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 61 const PetscInt **idx,*idx_i; 62 PetscInt *n,*w3,*w4,**data,len; 63 PetscErrorCode ierr; 64 PetscMPIInt size,rank,tag1,tag2,*w2,*w1,nrqr; 65 PetscInt Mbs,i,j,k,**rbuf,row,proc=-1,nrqs,msz,**outdat,**ptr; 66 PetscInt *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2,*d_p; 67 PetscMPIInt *onodes1,*olengths1,*onodes2,*olengths2; 68 PetscBT *table; 69 MPI_Comm comm; 70 MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 71 MPI_Status *s_status,*recv_status; 72 char *t_p; 73 74 PetscFunctionBegin; 75 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 76 size = c->size; 77 rank = c->rank; 78 Mbs = c->Mbs; 79 80 ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 81 ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 82 83 ierr = PetscMalloc2(imax+1,&idx,imax,&n);CHKERRQ(ierr); 84 85 for (i=0; i<imax; i++) { 86 ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr); 87 ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 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 ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); 93 for (i=0; i<imax; i++) { 94 ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* 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 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries"); 100 ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 101 w4[proc]++; 102 } 103 for (j=0; j<size; j++) { 104 if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 105 } 106 } 107 108 nrqs = 0; /* no of outgoing messages */ 109 msz = 0; /* total mesg length (for all proc */ 110 w1[rank] = 0; /* no mesg sent to itself */ 111 w3[rank] = 0; 112 for (i=0; i<size; i++) { 113 if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 114 } 115 /* pa - is list of processors to communicate with */ 116 ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); 117 for (i=0,j=0; i<size; i++) { 118 if (w1[i]) {pa[j] = i; j++;} 119 } 120 121 /* Each message would have a header = 1 + 2*(no of IS) + data */ 122 for (i=0; i<nrqs; i++) { 123 j = pa[i]; 124 w1[j] += w2[j] + 2*w3[j]; 125 msz += w1[j]; 126 } 127 128 /* Determine the number of messages to expect, their lengths, from from-ids */ 129 ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 130 ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 131 132 /* Now post the Irecvs corresponding to these messages */ 133 ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr); 134 135 /* Allocate Memory for outgoing messages */ 136 ierr = PetscMalloc4(size,&outdat,size,&ptr,msz,&tmp,size,&ctr);CHKERRQ(ierr); 137 ierr = PetscMemzero(outdat,size*sizeof(PetscInt*));CHKERRQ(ierr); 138 ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 139 { 140 PetscInt *iptr = tmp,ict = 0; 141 for (i=0; i<nrqs; i++) { 142 j = pa[i]; 143 iptr += ict; 144 outdat[j] = iptr; 145 ict = w1[j]; 146 } 147 } 148 149 /* Form the outgoing messages */ 150 /*plug in the headers*/ 151 for (i=0; i<nrqs; i++) { 152 j = pa[i]; 153 outdat[j][0] = 0; 154 ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 155 ptr[j] = outdat[j] + 2*w3[j] + 1; 156 } 157 158 /* Memory for doing local proc's work*/ 159 { 160 ierr = PetscCalloc5(imax,&table, imax,&data, imax,&isz, Mbs*imax,&d_p, (Mbs/PETSC_BITS_PER_BYTE+1)*imax,&t_p);CHKERRQ(ierr); 161 162 for (i=0; i<imax; i++) { 163 table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 164 data[i] = d_p + (Mbs)*i; 165 } 166 } 167 168 /* Parse the IS and update local tables and the outgoing buf with the data*/ 169 { 170 PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j; 171 PetscBT table_i; 172 173 for (i=0; i<imax; i++) { 174 ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 175 n_i = n[i]; 176 table_i = table[i]; 177 idx_i = idx[i]; 178 data_i = data[i]; 179 isz_i = isz[i]; 180 for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 181 row = idx_i[j]; 182 ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 183 if (proc != rank) { /* copy to the outgoing buffer */ 184 ctr[proc]++; 185 *ptr[proc] = row; 186 ptr[proc]++; 187 } else { /* Update the local table */ 188 if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 189 } 190 } 191 /* Update the headers for the current IS */ 192 for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 193 if ((ctr_j = ctr[j])) { 194 outdat_j = outdat[j]; 195 k = ++outdat_j[0]; 196 outdat_j[2*k] = ctr_j; 197 outdat_j[2*k-1] = i; 198 } 199 } 200 isz[i] = isz_i; 201 } 202 } 203 204 /* Now post the sends */ 205 ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 206 for (i=0; i<nrqs; ++i) { 207 j = pa[i]; 208 ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr); 209 } 210 211 /* No longer need the original indices*/ 212 for (i=0; i<imax; ++i) { 213 ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr); 214 } 215 ierr = PetscFree2(idx,n);CHKERRQ(ierr); 216 217 for (i=0; i<imax; ++i) { 218 ierr = ISDestroy(&is[i]);CHKERRQ(ierr); 219 } 220 221 /* Do Local work*/ 222 ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 223 224 /* Receive messages*/ 225 ierr = PetscMalloc1(nrqr+1,&recv_status);CHKERRQ(ierr); 226 if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);} 227 228 ierr = PetscMalloc1(nrqs+1,&s_status);CHKERRQ(ierr); 229 if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);} 230 231 /* Phase 1 sends are complete - deallocate buffers */ 232 ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr); 233 ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 234 235 ierr = PetscMalloc1(nrqr+1,&xdata);CHKERRQ(ierr); 236 ierr = PetscMalloc1(nrqr+1,&isz1);CHKERRQ(ierr); 237 ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 238 ierr = PetscFree(rbuf[0]);CHKERRQ(ierr); 239 ierr = PetscFree(rbuf);CHKERRQ(ierr); 240 241 /* Send the data back*/ 242 /* Do a global reduction to know the buffer space req for incoming messages*/ 243 { 244 PetscMPIInt *rw1; 245 246 ierr = PetscCalloc1(size,&rw1);CHKERRQ(ierr); 247 248 for (i=0; i<nrqr; ++i) { 249 proc = recv_status[i].MPI_SOURCE; 250 if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch"); 251 rw1[proc] = isz1[i]; 252 } 253 254 ierr = PetscFree(onodes1);CHKERRQ(ierr); 255 ierr = PetscFree(olengths1);CHKERRQ(ierr); 256 257 /* Determine the number of messages to expect, their lengths, from from-ids */ 258 ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 259 ierr = PetscFree(rw1);CHKERRQ(ierr); 260 } 261 /* Now post the Irecvs corresponding to these messages */ 262 ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 263 264 /* Now post the sends */ 265 ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 266 for (i=0; i<nrqr; ++i) { 267 j = recv_status[i].MPI_SOURCE; 268 ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr); 269 } 270 271 /* receive work done on other processors*/ 272 { 273 PetscMPIInt idex; 274 PetscInt is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 275 PetscBT table_i; 276 MPI_Status *status2; 277 278 ierr = PetscMalloc1(PetscMax(nrqr,nrqs)+1,&status2);CHKERRQ(ierr); 279 for (i=0; i<nrqs; ++i) { 280 ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr); 281 /* Process the message*/ 282 rbuf2_i = rbuf2[idex]; 283 ct1 = 2*rbuf2_i[0]+1; 284 jmax = rbuf2[idex][0]; 285 for (j=1; j<=jmax; j++) { 286 max = rbuf2_i[2*j]; 287 is_no = rbuf2_i[2*j-1]; 288 isz_i = isz[is_no]; 289 data_i = data[is_no]; 290 table_i = table[is_no]; 291 for (k=0; k<max; k++,ct1++) { 292 row = rbuf2_i[ct1]; 293 if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 294 } 295 isz[is_no] = isz_i; 296 } 297 } 298 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);} 299 ierr = PetscFree(status2);CHKERRQ(ierr); 300 } 301 302 for (i=0; i<imax; ++i) { 303 ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr); 304 } 305 306 307 ierr = PetscFree(onodes2);CHKERRQ(ierr); 308 ierr = PetscFree(olengths2);CHKERRQ(ierr); 309 310 ierr = PetscFree(pa);CHKERRQ(ierr); 311 ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr); 312 ierr = PetscFree(rbuf2);CHKERRQ(ierr); 313 ierr = PetscFree(s_waits1);CHKERRQ(ierr); 314 ierr = PetscFree(r_waits1);CHKERRQ(ierr); 315 ierr = PetscFree(s_waits2);CHKERRQ(ierr); 316 ierr = PetscFree(r_waits2);CHKERRQ(ierr); 317 ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr); 318 ierr = PetscFree(s_status);CHKERRQ(ierr); 319 ierr = PetscFree(recv_status);CHKERRQ(ierr); 320 ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 321 ierr = PetscFree(xdata);CHKERRQ(ierr); 322 ierr = PetscFree(isz1);CHKERRQ(ierr); 323 PetscFunctionReturn(0); 324 } 325 326 /* 327 MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 328 the work on the local processor. 329 330 Inputs: 331 C - MAT_MPIBAIJ; 332 imax - total no of index sets processed at a time; 333 table - an array of char - size = Mbs bits. 334 335 Output: 336 isz - array containing the count of the solution elements corresponding 337 to each index set; 338 data - pointer to the solutions 339 */ 340 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data) 341 { 342 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 343 Mat A = c->A,B = c->B; 344 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 345 PetscInt start,end,val,max,rstart,cstart,*ai,*aj; 346 PetscInt *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 347 PetscBT table_i; 348 349 PetscFunctionBegin; 350 rstart = c->rstartbs; 351 cstart = c->cstartbs; 352 ai = a->i; 353 aj = a->j; 354 bi = b->i; 355 bj = b->j; 356 garray = c->garray; 357 358 359 for (i=0; i<imax; i++) { 360 data_i = data[i]; 361 table_i = table[i]; 362 isz_i = isz[i]; 363 for (j=0,max=isz[i]; j<max; j++) { 364 row = data_i[j] - rstart; 365 start = ai[row]; 366 end = ai[row+1]; 367 for (k=start; k<end; k++) { /* Amat */ 368 val = aj[k] + cstart; 369 if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 370 } 371 start = bi[row]; 372 end = bi[row+1]; 373 for (k=start; k<end; k++) { /* Bmat */ 374 val = garray[bj[k]]; 375 if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 376 } 377 } 378 isz[i] = isz_i; 379 } 380 PetscFunctionReturn(0); 381 } 382 /* 383 MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 384 and return the output 385 386 Input: 387 C - the matrix 388 nrqr - no of messages being processed. 389 rbuf - an array of pointers to the recieved requests 390 391 Output: 392 xdata - array of messages to be sent back 393 isz1 - size of each message 394 395 For better efficiency perhaps we should malloc separately each xdata[i], 396 then if a remalloc is required we need only copy the data for that one row 397 rather than all previous rows as it is now where a single large chunck of 398 memory is used. 399 400 */ 401 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1) 402 { 403 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 404 Mat A = c->A,B = c->B; 405 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 406 PetscErrorCode ierr; 407 PetscInt rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 408 PetscInt row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 409 PetscInt val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr; 410 PetscInt *rbuf_i,kmax,rbuf_0; 411 PetscBT xtable; 412 413 PetscFunctionBegin; 414 Mbs = c->Mbs; 415 rstart = c->rstartbs; 416 cstart = c->cstartbs; 417 ai = a->i; 418 aj = a->j; 419 bi = b->i; 420 bj = b->j; 421 garray = c->garray; 422 423 424 for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 425 rbuf_i = rbuf[i]; 426 rbuf_0 = rbuf_i[0]; 427 ct += rbuf_0; 428 for (j=1; j<=rbuf_0; j++) total_sz += rbuf_i[2*j]; 429 } 430 431 if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs; 432 else max1 = 1; 433 mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 434 ierr = PetscMalloc1(mem_estimate,&xdata[0]);CHKERRQ(ierr); 435 ++no_malloc; 436 ierr = PetscBTCreate(Mbs,&xtable);CHKERRQ(ierr); 437 ierr = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr); 438 439 ct3 = 0; 440 for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 441 rbuf_i = rbuf[i]; 442 rbuf_0 = rbuf_i[0]; 443 ct1 = 2*rbuf_0+1; 444 ct2 = ct1; 445 ct3 += ct1; 446 for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 447 ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr); 448 oct2 = ct2; 449 kmax = rbuf_i[2*j]; 450 for (k=0; k<kmax; k++,ct1++) { 451 row = rbuf_i[ct1]; 452 if (!PetscBTLookupSet(xtable,row)) { 453 if (!(ct3 < mem_estimate)) { 454 new_estimate = (PetscInt)(1.5*mem_estimate)+1; 455 ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 456 ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 457 ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 458 xdata[0] = tmp; 459 mem_estimate = new_estimate; ++no_malloc; 460 for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 461 } 462 xdata[i][ct2++] = row; 463 ct3++; 464 } 465 } 466 for (k=oct2,max2=ct2; k<max2; k++) { 467 row = xdata[i][k] - rstart; 468 start = ai[row]; 469 end = ai[row+1]; 470 for (l=start; l<end; l++) { 471 val = aj[l] + cstart; 472 if (!PetscBTLookupSet(xtable,val)) { 473 if (!(ct3 < mem_estimate)) { 474 new_estimate = (PetscInt)(1.5*mem_estimate)+1; 475 ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 476 ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 477 ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 478 xdata[0] = tmp; 479 mem_estimate = new_estimate; ++no_malloc; 480 for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 481 } 482 xdata[i][ct2++] = val; 483 ct3++; 484 } 485 } 486 start = bi[row]; 487 end = bi[row+1]; 488 for (l=start; l<end; l++) { 489 val = garray[bj[l]]; 490 if (!PetscBTLookupSet(xtable,val)) { 491 if (!(ct3 < mem_estimate)) { 492 new_estimate = (PetscInt)(1.5*mem_estimate)+1; 493 ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 494 ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 495 ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 496 xdata[0] = tmp; 497 mem_estimate = new_estimate; ++no_malloc; 498 for (ctr =1; ctr <=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 499 } 500 xdata[i][ct2++] = val; 501 ct3++; 502 } 503 } 504 } 505 /* Update the header*/ 506 xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 507 xdata[i][2*j-1] = rbuf_i[2*j-1]; 508 } 509 xdata[i][0] = rbuf_0; 510 xdata[i+1] = xdata[i] + ct2; 511 isz1[i] = ct2; /* size of each message */ 512 } 513 ierr = PetscBTDestroy(&xtable);CHKERRQ(ierr); 514 ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr); 515 PetscFunctionReturn(0); 516 } 517 518 PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 519 { 520 IS *isrow_block,*iscol_block; 521 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 522 PetscErrorCode ierr; 523 PetscInt nmax,nstages_local,nstages,i,pos,max_no,ncol,nrow,N=C->cmap->N,bs=C->rmap->bs; 524 PetscBool colflag,*allcolumns,*allrows; 525 526 PetscFunctionBegin; 527 //printf("MatGetSubMatrices_MPIBAIJ scall %d, bs %d\n",scall,bs); 528 /* The compression and expansion should be avoided. Doesn't point 529 out errors, might change the indices, hence buggey */ 530 ierr = PetscMalloc2(ismax+1,&isrow_block,ismax+1,&iscol_block);CHKERRQ(ierr); 531 ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_block);CHKERRQ(ierr); 532 ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_block);CHKERRQ(ierr); 533 534 /* Check for special case: each processor gets entire matrix columns */ 535 ierr = PetscMalloc2(ismax+1,&allcolumns,ismax+1,&allrows);CHKERRQ(ierr); 536 for (i=0; i<ismax; i++) { 537 ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr); 538 ierr = ISGetLocalSize(iscol[i],&ncol);CHKERRQ(ierr); 539 if (colflag && ncol == C->cmap->N) allcolumns[i] = PETSC_TRUE; 540 else allcolumns[i] = PETSC_FALSE; 541 542 ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr); 543 ierr = ISGetLocalSize(isrow[i],&nrow);CHKERRQ(ierr); 544 if (colflag && nrow == C->rmap->N) allrows[i] = PETSC_TRUE; 545 else allrows[i] = PETSC_FALSE; 546 } 547 548 /* Allocate memory to hold all the submatrices */ 549 if (scall == MAT_INITIAL_MATRIX) { 550 ierr = PetscMalloc1(ismax+1,submat);CHKERRQ(ierr); 551 } 552 /* Determine the number of stages through which submatrices are done */ 553 nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 554 if (!nmax) nmax = 1; 555 nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0); 556 557 /* Make sure every processor loops through the nstages */ 558 ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr); 559 for (i=0,pos=0; i<nstages; i++) { 560 if (pos+nmax <= ismax) max_no = nmax; 561 else if (pos == ismax) max_no = 0; 562 else max_no = ismax-pos; 563 564 PetscBool isbaij; 565 ierr = PetscObjectTypeCompare((PetscObject)C,MATMPIBAIJ,&isbaij);CHKERRQ(ierr); 566 if (isbaij) { 567 ierr = MatGetSubMatrices_MPIBAIJ_local_new(C,max_no,isrow_block+pos,iscol_block+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr); 568 } else { 569 ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_block+pos,iscol_block+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr); 570 } 571 pos += max_no; 572 } 573 574 for (i=0; i<ismax; i++) { 575 ierr = ISDestroy(&isrow_block[i]);CHKERRQ(ierr); 576 ierr = ISDestroy(&iscol_block[i]);CHKERRQ(ierr); 577 } 578 ierr = PetscFree2(isrow_block,iscol_block);CHKERRQ(ierr); 579 ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr); 580 PetscFunctionReturn(0); 581 } 582 583 #if defined(PETSC_USE_CTABLE) 584 PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank) 585 { 586 PetscInt nGlobalNd = proc_gnode[size]; 587 PetscMPIInt fproc; 588 PetscErrorCode ierr; 589 590 PetscFunctionBegin; 591 ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr); 592 if (fproc > size) fproc = size; 593 while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) { 594 if (row < proc_gnode[fproc]) fproc--; 595 else fproc++; 596 } 597 *rank = fproc; 598 PetscFunctionReturn(0); 599 } 600 #endif 601 602 /* -------------------------------------------------------------------------*/ 603 /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 604 605 PetscErrorCode MatDestroy_MPIBAIJ_MatGetSubmatrices(Mat C) 606 { 607 PetscErrorCode ierr; 608 Mat_SeqBAIJ *c = (Mat_SeqBAIJ*)C->data; 609 Mat_SubMat *submatj = c->submatis1; 610 PetscInt i; 611 612 PetscFunctionBegin; 613 if (!submatj->id) { /* delete data that are linked only to submats[id=0] */ 614 ierr = PetscFree4(submatj->sbuf1,submatj->ptr,submatj->tmp,submatj->ctr);CHKERRQ(ierr); 615 616 for (i=0; i<submatj->nrqr; ++i) { 617 ierr = PetscFree(submatj->sbuf2[i]);CHKERRQ(ierr); 618 } 619 ierr = PetscFree3(submatj->sbuf2,submatj->req_size,submatj->req_source1);CHKERRQ(ierr); 620 621 if (submatj->rbuf1) { 622 ierr = PetscFree(submatj->rbuf1[0]);CHKERRQ(ierr); 623 ierr = PetscFree(submatj->rbuf1);CHKERRQ(ierr); 624 } 625 626 for (i=0; i<submatj->nrqs; ++i) { 627 ierr = PetscFree(submatj->rbuf3[i]);CHKERRQ(ierr); 628 } 629 ierr = PetscFree3(submatj->req_source2,submatj->rbuf2,submatj->rbuf3);CHKERRQ(ierr); 630 ierr = PetscFree(submatj->pa);CHKERRQ(ierr); 631 } 632 633 #if defined(PETSC_USE_CTABLE) 634 ierr = PetscTableDestroy((PetscTable*)&submatj->rmap);CHKERRQ(ierr); 635 if (submatj->cmap_loc) {ierr = PetscFree(submatj->cmap_loc);CHKERRQ(ierr);} 636 ierr = PetscFree(submatj->rmap_loc);CHKERRQ(ierr); 637 #else 638 ierr = PetscFree(submatj->rmap);CHKERRQ(ierr); 639 #endif 640 641 if (!submatj->allcolumns) { 642 #if defined(PETSC_USE_CTABLE) 643 ierr = PetscTableDestroy((PetscTable*)&submatj->cmap);CHKERRQ(ierr); 644 #else 645 ierr = PetscFree(submatj->cmap);CHKERRQ(ierr); 646 #endif 647 } 648 ierr = submatj->destroy(C);CHKERRQ(ierr); 649 ierr = PetscFree(submatj->row2proc);CHKERRQ(ierr); 650 651 ierr = PetscFree(submatj);CHKERRQ(ierr); 652 PetscFunctionReturn(0); 653 } 654 655 PetscErrorCode MatGetSubMatrices_MPIBAIJ_local_new(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats) 656 { 657 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 658 Mat A = c->A; 659 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*subc; 660 const PetscInt **icol,**irow; 661 PetscInt *nrow,*ncol,start; 662 PetscErrorCode ierr; 663 PetscMPIInt rank,size,tag0,tag2,tag3,tag4,*w1,*w2,*w3,*w4,nrqr; 664 PetscInt **sbuf1,**sbuf2,i,j,k,l,ct1,ct2,**rbuf1,row,proc=-1; 665 PetscInt nrqs=0,msz,**ptr=NULL,*req_size=NULL,*ctr=NULL,*pa,*tmp=NULL,tcol; 666 PetscInt **rbuf3=NULL,*req_source1=NULL,*req_source2,**sbuf_aj,**rbuf2=NULL,max1,max2; 667 PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax; 668 #if defined(PETSC_USE_CTABLE) 669 PetscTable *cmap,cmap_i=NULL,*rmap,rmap_i; 670 #else 671 PetscInt **cmap,*cmap_i=NULL,**rmap,*rmap_i; 672 #endif 673 const PetscInt *irow_i; 674 PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 675 MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 676 MPI_Request *r_waits4,*s_waits3,*s_waits4; 677 MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 678 MPI_Status *r_status3,*r_status4,*s_status4; 679 MPI_Comm comm; 680 PetscScalar **rbuf4,*rbuf4_i,**sbuf_aa,*vals,*mat_a,*imat_a,*sbuf_aa_i; 681 PetscMPIInt *onodes1,*olengths1,end; 682 PetscInt **row2proc,*row2proc_i,ilen_row,*imat_ilen,*imat_j,*imat_i; 683 Mat_SubMat **smats,*smat_i; 684 PetscBool *issorted,colflag,iscsorted=PETSC_TRUE; 685 PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,ilen; 686 687 PetscInt bs=C->rmap->bs,bs2=c->bs2,rstart = c->rstartbs; 688 PetscBool ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */ 689 PetscInt nzA,nzB,*a_i=a->i,*b_i=b->i,*a_j = a->j,*b_j = b->j,ctmp,imark,*cworkA,*cworkB; 690 PetscScalar *vworkA,*vworkB,*a_a = a->a,*b_a = b->a; 691 PetscInt cstart = c->cstartbs,*bmap = c->garray; 692 693 PetscFunctionBegin; 694 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 695 size = c->size; 696 rank = c->rank; 697 698 ierr = PetscMalloc5(ismax,&smats,ismax,&row2proc,ismax,&cmap,ismax,&rmap,ismax,&allcolumns);CHKERRQ(ierr); 699 ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,ismax,&issorted);CHKERRQ(ierr); 700 701 for (i=0; i<ismax; i++) { 702 ierr = ISSorted(iscol[i],&issorted[i]);CHKERRQ(ierr); 703 if (!issorted[i]) iscsorted = issorted[i]; 704 705 ierr = ISSorted(isrow[i],&issorted[i]);CHKERRQ(ierr); 706 707 ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 708 ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 709 710 /* Check for special case: allcolumn */ 711 ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr); 712 ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 713 if (colflag && ncol[i] == C->cmap->N) { 714 allcolumns[i] = PETSC_TRUE; 715 icol[i] = NULL; 716 } else { 717 allcolumns[i] = PETSC_FALSE; 718 ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 719 } 720 } 721 722 if (scall == MAT_REUSE_MATRIX) { 723 /* Assumes new rows are same length as the old rows */ 724 for (i=0; i<ismax; i++) { 725 subc = (Mat_SeqBAIJ*)(submats[i]->data); 726 if (subc->mbs != nrow[i] || subc->nbs != ncol[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 727 728 /* Initial matrix as if empty */ 729 ierr = PetscMemzero(subc->ilen,subc->mbs*sizeof(PetscInt));CHKERRQ(ierr); 730 731 /* Initial matrix as if empty */ 732 submats[i]->factortype = C->factortype; 733 734 smat_i = subc->submatis1; 735 smats[i] = smat_i; 736 737 nrqs = smat_i->nrqs; 738 nrqr = smat_i->nrqr; 739 rbuf1 = smat_i->rbuf1; 740 rbuf2 = smat_i->rbuf2; 741 rbuf3 = smat_i->rbuf3; 742 req_source2 = smat_i->req_source2; 743 744 sbuf1 = smat_i->sbuf1; 745 sbuf2 = smat_i->sbuf2; 746 ptr = smat_i->ptr; 747 tmp = smat_i->tmp; 748 ctr = smat_i->ctr; 749 750 pa = smat_i->pa; 751 req_size = smat_i->req_size; 752 req_source1 = smat_i->req_source1; 753 754 allcolumns[i] = smat_i->allcolumns; 755 row2proc[i] = smat_i->row2proc; 756 rmap[i] = smat_i->rmap; 757 cmap[i] = smat_i->cmap; 758 } 759 } else { /* scall == MAT_INITIAL_MATRIX */ 760 /* Get some new tags to keep the communication clean */ 761 ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 762 ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 763 764 /* evaluate communication - mesg to who, length of mesg, and buffer space 765 required. Based on this, buffers are allocated, and data copied into them*/ 766 ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); /* mesg size, initialize work vectors */ 767 768 for (i=0; i<ismax; i++) { 769 jmax = nrow[i]; 770 irow_i = irow[i]; 771 772 ierr = PetscMalloc1(jmax,&row2proc_i);CHKERRQ(ierr); 773 row2proc[i] = row2proc_i; 774 775 if (issorted[i]) proc = 0; 776 for (j=0; j<jmax; j++) { 777 if (!issorted[i]) proc = 0; 778 row = irow_i[j]; 779 while (row >= c->rangebs[proc+1]) proc++; 780 w4[proc]++; 781 row2proc_i[j] = proc; /* map row index to proc */ 782 } 783 for (j=0; j<size; j++) { 784 if (w4[j]) { w1[j] += w4[j]; w3[j]++; w4[j] = 0;} 785 } 786 } 787 788 nrqs = 0; /* no of outgoing messages */ 789 msz = 0; /* total mesg length (for all procs) */ 790 w1[rank] = 0; /* no mesg sent to self */ 791 w3[rank] = 0; 792 for (i=0; i<size; i++) { 793 if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 794 } 795 ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/ 796 for (i=0,j=0; i<size; i++) { 797 if (w1[i]) { pa[j] = i; j++; } 798 } 799 800 /* Each message would have a header = 1 + 2*(no of IS) + data */ 801 for (i=0; i<nrqs; i++) { 802 j = pa[i]; 803 w1[j] += w2[j] + 2* w3[j]; 804 msz += w1[j]; 805 } 806 ierr = PetscInfo2(0,"Number of outgoing messages %D Total message length %D\n",nrqs,msz);CHKERRQ(ierr); 807 808 /* Determine the number of messages to expect, their lengths, from from-ids */ 809 ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 810 ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 811 812 /* Now post the Irecvs corresponding to these messages */ 813 tag0 = ((PetscObject)C)->tag; 814 ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 815 816 ierr = PetscFree(onodes1);CHKERRQ(ierr); 817 ierr = PetscFree(olengths1);CHKERRQ(ierr); 818 819 /* Allocate Memory for outgoing messages */ 820 ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr); 821 ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr); 822 ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 823 824 { 825 PetscInt *iptr = tmp; 826 k = 0; 827 for (i=0; i<nrqs; i++) { 828 j = pa[i]; 829 iptr += k; 830 sbuf1[j] = iptr; 831 k = w1[j]; 832 } 833 } 834 835 /* Form the outgoing messages. Initialize the header space */ 836 for (i=0; i<nrqs; i++) { 837 j = pa[i]; 838 sbuf1[j][0] = 0; 839 ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 840 ptr[j] = sbuf1[j] + 2*w3[j] + 1; 841 } 842 843 /* Parse the isrow and copy data into outbuf */ 844 for (i=0; i<ismax; i++) { 845 row2proc_i = row2proc[i]; 846 ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 847 irow_i = irow[i]; 848 jmax = nrow[i]; 849 for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 850 proc = row2proc_i[j]; 851 if (proc != rank) { /* copy to the outgoing buf*/ 852 ctr[proc]++; 853 *ptr[proc] = irow_i[j]; 854 ptr[proc]++; 855 } 856 } 857 /* Update the headers for the current IS */ 858 for (j=0; j<size; j++) { /* Can Optimise this loop too */ 859 if ((ctr_j = ctr[j])) { 860 sbuf1_j = sbuf1[j]; 861 k = ++sbuf1_j[0]; 862 sbuf1_j[2*k] = ctr_j; 863 sbuf1_j[2*k-1] = i; 864 } 865 } 866 } 867 868 /* Now post the sends */ 869 ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 870 for (i=0; i<nrqs; ++i) { 871 j = pa[i]; 872 ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 873 } 874 875 /* Post Receives to capture the buffer size */ 876 ierr = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr); 877 ierr = PetscMalloc3(nrqs+1,&req_source2,nrqs+1,&rbuf2,nrqs+1,&rbuf3);CHKERRQ(ierr); 878 rbuf2[0] = tmp + msz; 879 for (i=1; i<nrqs; ++i) { 880 rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 881 } 882 for (i=0; i<nrqs; ++i) { 883 j = pa[i]; 884 ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag2,comm,r_waits2+i);CHKERRQ(ierr); 885 } 886 887 /* Send to other procs the buf size they should allocate */ 888 /* Receive messages*/ 889 ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 890 ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr); 891 ierr = PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);CHKERRQ(ierr); 892 { 893 PetscInt *sAi = a->i,*sBi = b->i,id; 894 PetscInt *sbuf2_i; 895 896 ierr = MPI_Waitall(nrqr,r_waits1,r_status1);CHKERRQ(ierr); 897 for (i=0; i<nrqr; ++i) { 898 req_size[i] = 0; 899 rbuf1_i = rbuf1[i]; 900 start = 2*rbuf1_i[0] + 1; 901 ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 902 ierr = PetscMalloc1(end+1,&sbuf2[i]);CHKERRQ(ierr); 903 sbuf2_i = sbuf2[i]; 904 for (j=start; j<end; j++) { 905 id = rbuf1_i[j] - rstart; 906 ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 907 sbuf2_i[j] = ncols; 908 req_size[i] += ncols; 909 } 910 req_source1[i] = r_status1[i].MPI_SOURCE; 911 /* form the header */ 912 sbuf2_i[0] = req_size[i]; 913 for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j]; 914 915 ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);CHKERRQ(ierr); 916 } 917 } 918 ierr = PetscFree(r_status1);CHKERRQ(ierr); 919 ierr = PetscFree(r_waits1);CHKERRQ(ierr); 920 ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 921 922 /* Receive messages*/ 923 ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr); 924 ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr); 925 926 ierr = MPI_Waitall(nrqs,r_waits2,r_status2);CHKERRQ(ierr); 927 for (i=0; i<nrqs; ++i) { 928 ierr = PetscMalloc1(rbuf2[i][0]+1,&rbuf3[i]);CHKERRQ(ierr); 929 req_source2[i] = r_status2[i].MPI_SOURCE; 930 ierr = MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);CHKERRQ(ierr); 931 } 932 ierr = PetscFree(r_status2);CHKERRQ(ierr); 933 ierr = PetscFree(r_waits2);CHKERRQ(ierr); 934 935 /* Wait on sends1 and sends2 */ 936 ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr); 937 ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr); 938 939 if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);} 940 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);} 941 ierr = PetscFree(s_status1);CHKERRQ(ierr); 942 ierr = PetscFree(s_status2);CHKERRQ(ierr); 943 ierr = PetscFree(s_waits1);CHKERRQ(ierr); 944 ierr = PetscFree(s_waits2);CHKERRQ(ierr); 945 946 /* Now allocate sending buffers for a->j, and send them off */ 947 ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr); 948 for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 949 ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr); 950 for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 951 952 ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr); 953 { 954 955 for (i=0; i<nrqr; i++) { 956 rbuf1_i = rbuf1[i]; 957 sbuf_aj_i = sbuf_aj[i]; 958 ct1 = 2*rbuf1_i[0] + 1; 959 ct2 = 0; 960 for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 961 kmax = rbuf1[i][2*j]; 962 for (k=0; k<kmax; k++,ct1++) { 963 row = rbuf1_i[ct1] - rstart; 964 nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 965 ncols = nzA + nzB; 966 cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 967 968 /* load the column indices for this row into cols */ 969 cols = sbuf_aj_i + ct2; 970 for (l=0; l<nzB; l++) { 971 if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 972 else break; 973 } 974 imark = l; 975 for (l=0; l<nzA; l++) {cols[imark+l] = cstart + cworkA[l];} 976 for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 977 ct2 += ncols; 978 } 979 } 980 ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);CHKERRQ(ierr); 981 } 982 } 983 ierr = PetscMalloc2(nrqs+1,&r_status3,nrqr+1,&s_status3);CHKERRQ(ierr); 984 985 /* create col map: global col of C -> local col of submatrices */ 986 const PetscInt *icol_i; 987 #if defined(PETSC_USE_CTABLE) 988 for (i=0; i<ismax; i++) { 989 if (!allcolumns[i]) { 990 ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr); 991 992 jmax = ncol[i]; 993 icol_i = icol[i]; 994 cmap_i = cmap[i]; 995 for (j=0; j<jmax; j++) { 996 ierr = PetscTableAdd(cmap[i],icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 997 } 998 } else cmap[i] = NULL; 999 } 1000 #else 1001 for (i=0; i<ismax; i++) { 1002 if (!allcolumns[i]) { 1003 ierr = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr); 1004 jmax = ncol[i]; 1005 icol_i = icol[i]; 1006 cmap_i = cmap[i]; 1007 for (j=0; j<jmax; j++) cmap_i[icol_i[j]] = j+1; 1008 } else cmap[i] = NULL; 1009 } 1010 #endif 1011 1012 /* Create lens which is required for MatCreate... */ 1013 for (i=0,j=0; i<ismax; i++) j += nrow[i]; 1014 ierr = PetscMalloc1(ismax,&lens);CHKERRQ(ierr); 1015 1016 if (ismax) { 1017 ierr = PetscCalloc1(j,&lens[0]);CHKERRQ(ierr); 1018 } 1019 for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1]; 1020 1021 /* Update lens from local data */ 1022 for (i=0; i<ismax; i++) { 1023 row2proc_i = row2proc[i]; 1024 jmax = nrow[i]; 1025 if (!allcolumns[i]) cmap_i = cmap[i]; 1026 irow_i = irow[i]; 1027 lens_i = lens[i]; 1028 for (j=0; j<jmax; j++) { 1029 row = irow_i[j]; /* global blocked row of C */ 1030 proc = row2proc_i[j]; 1031 if (proc == rank) { 1032 /* Get indices from matA and then from matB */ 1033 #if defined(PETSC_USE_CTABLE) 1034 PetscInt tt; 1035 #endif 1036 row = row - rstart; 1037 nzA = a_i[row+1] - a_i[row]; 1038 nzB = b_i[row+1] - b_i[row]; 1039 cworkA = a_j + a_i[row]; 1040 cworkB = b_j + b_i[row]; 1041 1042 if (!allcolumns[i]) { 1043 #if defined(PETSC_USE_CTABLE) 1044 for (k=0; k<nzA; k++) { 1045 ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1046 if (tt) lens_i[j]++; 1047 } 1048 for (k=0; k<nzB; k++) { 1049 ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1050 if (tt) lens_i[j]++; 1051 } 1052 1053 #else 1054 for (k=0; k<nzA; k++) { 1055 if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1056 } 1057 for (k=0; k<nzB; k++) { 1058 if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1059 } 1060 #endif 1061 } else { /* allcolumns */ 1062 lens_i[j] = nzA + nzB; 1063 } 1064 } 1065 } 1066 } 1067 1068 /* Create row map: global row of C -> local row of submatrices */ 1069 #if defined(PETSC_USE_CTABLE) 1070 for (i=0; i<ismax; i++) { 1071 ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr); 1072 irow_i = irow[i]; 1073 jmax = nrow[i]; 1074 for (j=0; j<jmax; j++) { 1075 ierr = PetscTableAdd(rmap[i],irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1076 } 1077 } 1078 #else 1079 for (i=0; i<ismax; i++) { 1080 ierr = PetscCalloc1(c->Mbs,&rmap[i]);CHKERRQ(ierr); 1081 rmap_i = rmap[i]; 1082 irow_i = irow[i]; 1083 jmax = nrow[i]; 1084 for (j=0; j<jmax; j++) { 1085 rmap_i[irow_i[j]] = j; 1086 } 1087 } 1088 #endif 1089 1090 /* Update lens from offproc data */ 1091 { 1092 PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 1093 1094 ierr = MPI_Waitall(nrqs,r_waits3,r_status3);CHKERRQ(ierr); 1095 for (tmp2=0; tmp2<nrqs; tmp2++) { 1096 sbuf1_i = sbuf1[pa[tmp2]]; 1097 jmax = sbuf1_i[0]; 1098 ct1 = 2*jmax+1; 1099 ct2 = 0; 1100 rbuf2_i = rbuf2[tmp2]; 1101 rbuf3_i = rbuf3[tmp2]; 1102 for (j=1; j<=jmax; j++) { 1103 is_no = sbuf1_i[2*j-1]; 1104 max1 = sbuf1_i[2*j]; 1105 lens_i = lens[is_no]; 1106 if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1107 rmap_i = rmap[is_no]; 1108 for (k=0; k<max1; k++,ct1++) { 1109 #if defined(PETSC_USE_CTABLE) 1110 ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 1111 row--; 1112 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1113 #else 1114 row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1115 #endif 1116 max2 = rbuf2_i[ct1]; 1117 for (l=0; l<max2; l++,ct2++) { 1118 if (!allcolumns[is_no]) { 1119 #if defined(PETSC_USE_CTABLE) 1120 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1121 #else 1122 tcol = cmap_i[rbuf3_i[ct2]]; 1123 #endif 1124 if (tcol) lens_i[row]++; 1125 } else { /* allcolumns */ 1126 lens_i[row]++; /* lens_i[row] += max2 ? */ 1127 } 1128 } 1129 } 1130 } 1131 } 1132 } 1133 ierr = PetscFree(r_waits3);CHKERRQ(ierr); 1134 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);} 1135 ierr = PetscFree2(r_status3,s_status3);CHKERRQ(ierr); 1136 ierr = PetscFree(s_waits3);CHKERRQ(ierr); 1137 1138 /* Create the submatrices */ 1139 for (i=0; i<ismax; i++) { 1140 PetscInt bs_tmp; 1141 if (ijonly) bs_tmp = 1; 1142 else bs_tmp = bs; 1143 1144 ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 1145 ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1146 1147 ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 1148 ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); 1149 ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */ 1150 1151 /* create struct Mat_SubMat and attached it to submat */ 1152 ierr = PetscNew(&smat_i);CHKERRQ(ierr); 1153 subc = (Mat_SeqBAIJ*)submats[i]->data; 1154 subc->submatis1 = smat_i; 1155 smats[i] = smat_i; 1156 1157 smat_i->destroy = submats[i]->ops->destroy; 1158 submats[i]->ops->destroy = MatDestroy_MPIBAIJ_MatGetSubmatrices; 1159 submats[i]->factortype = C->factortype; 1160 1161 smat_i->id = i; 1162 smat_i->nrqs = nrqs; 1163 smat_i->nrqr = nrqr; 1164 smat_i->rbuf1 = rbuf1; 1165 smat_i->rbuf2 = rbuf2; 1166 smat_i->rbuf3 = rbuf3; 1167 smat_i->sbuf2 = sbuf2; 1168 smat_i->req_source2 = req_source2; 1169 1170 smat_i->sbuf1 = sbuf1; 1171 smat_i->ptr = ptr; 1172 smat_i->tmp = tmp; 1173 smat_i->ctr = ctr; 1174 1175 smat_i->pa = pa; 1176 smat_i->req_size = req_size; 1177 smat_i->req_source1 = req_source1; 1178 1179 smat_i->allcolumns = allcolumns[i]; 1180 smat_i->singleis = PETSC_FALSE; 1181 smat_i->row2proc = row2proc[i]; 1182 smat_i->rmap = rmap[i]; 1183 smat_i->cmap = cmap[i]; 1184 } 1185 1186 if (ismax) {ierr = PetscFree(lens[0]);CHKERRQ(ierr);} 1187 ierr = PetscFree(lens);CHKERRQ(ierr); 1188 ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 1189 ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 1190 1191 } /* endof scall == MAT_INITIAL_MATRIX */ 1192 1193 /* Post recv matrix values */ 1194 ierr = PetscObjectGetNewTag((PetscObject)C,&tag4);CHKERRQ(ierr); 1195 ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr); 1196 ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr); 1197 ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr); 1198 ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr); 1199 for (i=0; i<nrqs; ++i) { 1200 ierr = PetscMalloc1(rbuf2[i][0]*bs2,&rbuf4[i]);CHKERRQ(ierr); 1201 ierr = MPI_Irecv(rbuf4[i],rbuf2[i][0]*bs2,MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);CHKERRQ(ierr); 1202 } 1203 1204 /* Allocate sending buffers for a->a, and send them off */ 1205 ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr); 1206 for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1207 1208 ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr); 1209 for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 1210 1211 ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr); 1212 1213 for (i=0; i<nrqr; i++) { 1214 rbuf1_i = rbuf1[i]; 1215 sbuf_aa_i = sbuf_aa[i]; 1216 ct1 = 2*rbuf1_i[0]+1; 1217 ct2 = 0; 1218 for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1219 kmax = rbuf1_i[2*j]; 1220 for (k=0; k<kmax; k++,ct1++) { 1221 row = rbuf1_i[ct1] - rstart; 1222 nzA = a_i[row+1] - a_i[row]; 1223 nzB = b_i[row+1] - b_i[row]; 1224 ncols = nzA + nzB; 1225 cworkB = b_j + b_i[row]; 1226 vworkA = a_a + a_i[row]*bs2; 1227 vworkB = b_a + b_i[row]*bs2; 1228 1229 /* load the column values for this row into vals*/ 1230 vals = sbuf_aa_i+ct2*bs2; 1231 for (l=0; l<nzB; l++) { 1232 if ((bmap[cworkB[l]]) < cstart) { 1233 ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1234 } else break; 1235 } 1236 imark = l; 1237 for (l=0; l<nzA; l++) { 1238 ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);//error in proc[1]??? 1239 } 1240 for (l=imark; l<nzB; l++) { 1241 ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1242 } 1243 1244 ct2 += ncols; 1245 } 1246 } 1247 ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);CHKERRQ(ierr); 1248 } 1249 1250 if (!ismax) { 1251 ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr); 1252 ierr = PetscFree(rbuf1);CHKERRQ(ierr); 1253 } 1254 1255 /* Assemble the matrices */ 1256 /* First assemble the local rows */ 1257 for (i=0; i<ismax; i++) { 1258 row2proc_i = row2proc[i]; 1259 subc = (Mat_SeqBAIJ*)submats[i]->data; 1260 imat_ilen = subc->ilen; 1261 imat_j = subc->j; 1262 imat_i = subc->i; 1263 imat_a = subc->a; 1264 1265 if (!allcolumns[i]) cmap_i = cmap[i]; 1266 rmap_i = rmap[i]; 1267 irow_i = irow[i]; 1268 jmax = nrow[i]; 1269 for (j=0; j<jmax; j++) { 1270 row = irow_i[j]; 1271 proc = row2proc_i[j]; 1272 1273 if (proc == rank) { 1274 1275 row = row - rstart; 1276 nzA = a_i[row+1] - a_i[row]; 1277 nzB = b_i[row+1] - b_i[row]; 1278 cworkA = a_j + a_i[row]; 1279 cworkB = b_j + b_i[row]; 1280 if (!ijonly) { 1281 vworkA = a_a + a_i[row]*bs2; 1282 vworkB = b_a + b_i[row]*bs2; 1283 } 1284 #if defined(PETSC_USE_CTABLE) 1285 ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr); 1286 row--; 1287 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1288 #else 1289 row = rmap_i[row + rstart]; 1290 #endif 1291 mat_i = imat_i[row]; 1292 if (!ijonly) mat_a = imat_a + mat_i*bs2; 1293 mat_j = imat_j + mat_i; 1294 ilen_row = imat_ilen[row]; 1295 1296 /* load the column indices for this row into cols*/ 1297 if (!allcolumns[i]) { 1298 for (l=0; l<nzB; l++) { 1299 if ((ctmp = bmap[cworkB[l]]) < cstart) { 1300 #if defined(PETSC_USE_CTABLE) 1301 ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr); 1302 if (tcol) { 1303 #else 1304 if ((tcol = cmap_i[ctmp])) { 1305 #endif 1306 *mat_j++ = tcol - 1; 1307 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1308 mat_a += bs2; 1309 ilen_row++; 1310 } 1311 } else break; 1312 } 1313 imark = l; 1314 for (l=0; l<nzA; l++) { 1315 #if defined(PETSC_USE_CTABLE) 1316 ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1317 if (tcol) { 1318 #else 1319 if ((tcol = cmap_i[cstart + cworkA[l]])) { 1320 #endif 1321 *mat_j++ = tcol - 1; 1322 if (!ijonly) { 1323 ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1324 mat_a += bs2; 1325 } 1326 ilen_row++; 1327 } 1328 } 1329 for (l=imark; l<nzB; l++) { 1330 #if defined(PETSC_USE_CTABLE) 1331 ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1332 if (tcol) { 1333 #else 1334 if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1335 #endif 1336 *mat_j++ = tcol - 1; 1337 if (!ijonly) { 1338 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1339 mat_a += bs2; 1340 } 1341 ilen_row++; 1342 } 1343 } 1344 } else { /* allcolumns */ 1345 for (l=0; l<nzB; l++) { 1346 if ((ctmp = bmap[cworkB[l]]) < cstart) { 1347 *mat_j++ = ctmp; 1348 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1349 mat_a += bs2; 1350 ilen_row++; 1351 } else break; 1352 } 1353 imark = l; 1354 for (l=0; l<nzA; l++) { 1355 *mat_j++ = cstart+cworkA[l]; 1356 if (!ijonly) { 1357 ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1358 mat_a += bs2; 1359 } 1360 ilen_row++; 1361 } 1362 for (l=imark; l<nzB; l++) { 1363 *mat_j++ = bmap[cworkB[l]]; 1364 if (!ijonly) { 1365 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1366 mat_a += bs2; 1367 } 1368 ilen_row++; 1369 } 1370 } 1371 imat_ilen[row] = ilen_row; 1372 } 1373 } 1374 } 1375 1376 /* Now assemble the off proc rows */ 1377 ierr = MPI_Waitall(nrqs,r_waits4,r_status4);CHKERRQ(ierr); 1378 for (tmp2=0; tmp2<nrqs; tmp2++) { 1379 sbuf1_i = sbuf1[pa[tmp2]]; 1380 jmax = sbuf1_i[0]; 1381 ct1 = 2*jmax + 1; 1382 ct2 = 0; 1383 rbuf2_i = rbuf2[tmp2]; 1384 rbuf3_i = rbuf3[tmp2]; 1385 rbuf4_i = rbuf4[tmp2]; 1386 for (j=1; j<=jmax; j++) { 1387 is_no = sbuf1_i[2*j-1]; 1388 rmap_i = rmap[is_no]; 1389 if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1390 subc = (Mat_SeqBAIJ*)submats[is_no]->data; 1391 imat_ilen = subc->ilen; 1392 imat_j = subc->j; 1393 imat_i = subc->i; 1394 imat_a = subc->a; 1395 max1 = sbuf1_i[2*j]; 1396 for (k=0; k<max1; k++,ct1++) { 1397 row = sbuf1_i[ct1]; 1398 #if defined(PETSC_USE_CTABLE) 1399 ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr); 1400 row--; 1401 #else 1402 row = rmap_i[row]; 1403 #endif 1404 ilen = imat_ilen[row]; 1405 mat_i = imat_i[row]; 1406 mat_a = imat_a + mat_i; 1407 mat_j = imat_j + mat_i; 1408 max2 = rbuf2_i[ct1]; 1409 if (!allcolumns[is_no]) { 1410 for (l=0; l<max2; l++,ct2++) { 1411 #if defined(PETSC_USE_CTABLE) 1412 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1413 #else 1414 tcol = cmap_i[rbuf3_i[ct2]]; 1415 #endif 1416 if (tcol) { 1417 *mat_j++ = tcol - 1; 1418 *mat_a++ = rbuf4_i[ct2]; 1419 ilen++; 1420 } 1421 } 1422 } else { /* allcolumns */ 1423 for (l=0; l<max2; l++,ct2++) { 1424 *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */ 1425 *mat_a++ = rbuf4_i[ct2]; 1426 ilen++; 1427 } 1428 } 1429 imat_ilen[row] = ilen; 1430 } 1431 } 1432 } 1433 1434 if (!iscsorted) { /* sort column indices of the rows */ 1435 for (i=0; i<ismax; i++) { 1436 subc = (Mat_SeqBAIJ*)submats[i]->data; 1437 imat_j = subc->j; 1438 imat_i = subc->i; 1439 imat_a = subc->a; 1440 imat_ilen = subc->ilen; 1441 1442 if (allcolumns[i]) continue; 1443 jmax = nrow[i]; 1444 for (j=0; j<jmax; j++) { 1445 PetscInt ilen; 1446 1447 mat_i = imat_i[j]; 1448 mat_a = imat_a + mat_i; 1449 mat_j = imat_j + mat_i; 1450 ilen = imat_ilen[j]; 1451 ierr = PetscSortIntWithScalarArray(ilen,mat_j,mat_a);CHKERRQ(ierr); 1452 } 1453 } 1454 } 1455 1456 ierr = PetscFree(r_status4);CHKERRQ(ierr); 1457 ierr = PetscFree(r_waits4);CHKERRQ(ierr); 1458 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);} 1459 ierr = PetscFree(s_waits4);CHKERRQ(ierr); 1460 ierr = PetscFree(s_status4);CHKERRQ(ierr); 1461 1462 /* Restore the indices */ 1463 for (i=0; i<ismax; i++) { 1464 ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 1465 if (!allcolumns[i]) { 1466 ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 1467 } 1468 } 1469 1470 for (i=0; i<ismax; i++) { 1471 ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1472 ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1473 } 1474 1475 /* Destroy allocated memory */ 1476 if (!ismax) { 1477 ierr = PetscFree(pa);CHKERRQ(ierr); 1478 1479 ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr); 1480 for (i=0; i<nrqr; ++i) { 1481 ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 1482 } 1483 for (i=0; i<nrqs; ++i) { 1484 ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 1485 } 1486 1487 ierr = PetscFree3(sbuf2,req_size,req_source1);CHKERRQ(ierr); 1488 ierr = PetscFree3(req_source2,rbuf2,rbuf3);CHKERRQ(ierr); 1489 } 1490 1491 ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 1492 ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 1493 ierr = PetscFree5(irow,icol,nrow,ncol,issorted);CHKERRQ(ierr); 1494 1495 for (i=0; i<nrqs; ++i) { 1496 ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr); 1497 } 1498 ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1499 1500 ierr = PetscFree5(smats,row2proc,cmap,rmap,allcolumns);CHKERRQ(ierr); 1501 PetscFunctionReturn(0); 1502 } 1503 1504 //------------ endof new ------------- 1505 1506 PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats) 1507 { 1508 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 1509 Mat A = c->A; 1510 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat; 1511 const PetscInt **irow,**icol,*irow_i; 1512 PetscInt *nrow,*ncol,*w3,*w4,start; 1513 PetscErrorCode ierr; 1514 PetscMPIInt size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc; 1515 PetscInt **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row; 1516 PetscInt nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol; 1517 PetscInt **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2; 1518 PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax; 1519 PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 1520 PetscInt bs =C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB; 1521 PetscInt cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 1522 PetscInt *bmap = c->garray,ctmp,rstart=c->rstartbs; 1523 MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3,*s_waits3; 1524 MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2,*r_status3; 1525 MPI_Comm comm; 1526 PetscBool flag; 1527 PetscMPIInt *onodes1,*olengths1; 1528 PetscBool ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */ 1529 1530 /* variables below are used for the matrix numerical values - case of !ijonly */ 1531 MPI_Request *r_waits4,*s_waits4; 1532 MPI_Status *r_status4,*s_status4; 1533 MatScalar **rbuf4,**sbuf_aa,*vals,*mat_a = NULL,*sbuf_aa_i,*vworkA = NULL,*vworkB = NULL; 1534 MatScalar *a_a=a->a,*b_a=b->a; 1535 1536 #if defined(PETSC_USE_CTABLE) 1537 PetscInt tt; 1538 PetscTable *rmap,*cmap,rmap_i,cmap_i=NULL; 1539 #else 1540 PetscInt **cmap,*cmap_i=NULL,*rtable,*rmap_i,**rmap, Mbs = c->Mbs; 1541 #endif 1542 1543 PetscFunctionBegin; 1544 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1545 tag0 = ((PetscObject)C)->tag; 1546 size = c->size; 1547 rank = c->rank; 1548 //if (!rank) printf("MatGetSubMatrices_MPIBAIJ scall %d, bs %d\n",scall,bs); 1549 1550 /* Get some new tags to keep the communication clean */ 1551 ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 1552 ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 1553 ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 1554 1555 #if defined(PETSC_USE_CTABLE) 1556 ierr = PetscMalloc4(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol);CHKERRQ(ierr); 1557 #else 1558 ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,Mbs+1,&rtable);CHKERRQ(ierr); 1559 /* Create hash table for the mapping :row -> proc*/ 1560 for (i=0,j=0; i<size; i++) { 1561 jmax = C->rmap->range[i+1]/bs; 1562 for (; j<jmax; j++) rtable[j] = i; 1563 } 1564 #endif 1565 1566 for (i=0; i<ismax; i++) { 1567 if (allrows[i]) { 1568 irow[i] = NULL; 1569 nrow[i] = C->rmap->N/bs; 1570 } else { 1571 ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 1572 ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 1573 } 1574 1575 if (allcolumns[i]) { 1576 icol[i] = NULL; 1577 ncol[i] = C->cmap->N/bs; 1578 } else { 1579 ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 1580 ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 1581 } 1582 } 1583 1584 /* evaluate communication - mesg to who,length of mesg,and buffer space 1585 required. Based on this, buffers are allocated, and data copied into them*/ 1586 ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); 1587 for (i=0; i<ismax; i++) { 1588 ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 1589 jmax = nrow[i]; 1590 irow_i = irow[i]; 1591 for (j=0; j<jmax; j++) { 1592 if (allrows[i]) row = j; 1593 else row = irow_i[j]; 1594 1595 #if defined(PETSC_USE_CTABLE) 1596 ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1597 #else 1598 proc = rtable[row]; 1599 #endif 1600 w4[proc]++; 1601 } 1602 for (j=0; j<size; j++) { 1603 if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 1604 } 1605 } 1606 1607 nrqs = 0; /* no of outgoing messages */ 1608 msz = 0; /* total mesg length for all proc */ 1609 w1[rank] = 0; /* no mesg sent to intself */ 1610 w3[rank] = 0; 1611 for (i=0; i<size; i++) { 1612 if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 1613 } 1614 ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/ 1615 for (i=0,j=0; i<size; i++) { 1616 if (w1[i]) { pa[j] = i; j++; } 1617 } 1618 1619 /* Each message would have a header = 1 + 2*(no of IS) + data */ 1620 for (i=0; i<nrqs; i++) { 1621 j = pa[i]; 1622 w1[j] += w2[j] + 2* w3[j]; 1623 msz += w1[j]; 1624 } 1625 1626 /* Determine the number of messages to expect, their lengths, from from-ids */ 1627 ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 1628 ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 1629 1630 /* Now post the Irecvs corresponding to these messages */ 1631 ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 1632 1633 ierr = PetscFree(onodes1);CHKERRQ(ierr); 1634 ierr = PetscFree(olengths1);CHKERRQ(ierr); 1635 1636 /* Allocate Memory for outgoing messages */ 1637 ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr); 1638 ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr); 1639 ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 1640 { 1641 PetscInt *iptr = tmp,ict = 0; 1642 for (i=0; i<nrqs; i++) { 1643 j = pa[i]; 1644 iptr += ict; 1645 sbuf1[j] = iptr; 1646 ict = w1[j]; 1647 } 1648 } 1649 1650 /* Form the outgoing messages */ 1651 /* Initialise the header space */ 1652 for (i=0; i<nrqs; i++) { 1653 j = pa[i]; 1654 sbuf1[j][0] = 0; 1655 ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 1656 ptr[j] = sbuf1[j] + 2*w3[j] + 1; 1657 } 1658 1659 /* Parse the isrow and copy data into outbuf */ 1660 for (i=0; i<ismax; i++) { 1661 ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 1662 irow_i = irow[i]; 1663 jmax = nrow[i]; 1664 for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 1665 if (allrows[i]) row = j; 1666 else row = irow_i[j]; 1667 1668 #if defined(PETSC_USE_CTABLE) 1669 ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1670 #else 1671 proc = rtable[row]; 1672 #endif 1673 if (proc != rank) { /* copy to the outgoing buf*/ 1674 ctr[proc]++; 1675 *ptr[proc] = row; 1676 ptr[proc]++; 1677 } 1678 } 1679 /* Update the headers for the current IS */ 1680 for (j=0; j<size; j++) { /* Can Optimise this loop too */ 1681 if ((ctr_j = ctr[j])) { 1682 sbuf1_j = sbuf1[j]; 1683 k = ++sbuf1_j[0]; 1684 sbuf1_j[2*k] = ctr_j; 1685 sbuf1_j[2*k-1] = i; 1686 } 1687 } 1688 } 1689 1690 /* Now post the sends */ 1691 ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 1692 for (i=0; i<nrqs; ++i) { 1693 j = pa[i]; 1694 ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 1695 } 1696 1697 /* Post Recieves to capture the buffer size */ 1698 ierr = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr); 1699 ierr = PetscMalloc1(nrqs+1,&rbuf2);CHKERRQ(ierr); 1700 rbuf2[0] = tmp + msz; 1701 for (i=1; i<nrqs; ++i) { 1702 rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 1703 } 1704 for (i=0; i<nrqs; ++i) { 1705 j = pa[i]; 1706 ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 1707 } 1708 1709 /* Send to other procs the buf size they should allocate */ 1710 1711 /* Receive messages*/ 1712 ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 1713 ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr); 1714 ierr = PetscMalloc3(nrqr+1,&sbuf2,nrqr,&req_size,nrqr,&req_source);CHKERRQ(ierr); 1715 { 1716 Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data; 1717 PetscInt *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i; 1718 1719 for (i=0; i<nrqr; ++i) { 1720 ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr); 1721 1722 req_size[idex] = 0; 1723 rbuf1_i = rbuf1[idex]; 1724 start = 2*rbuf1_i[0] + 1; 1725 ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 1726 ierr = PetscMalloc1(end,&sbuf2[idex]);CHKERRQ(ierr); 1727 sbuf2_i = sbuf2[idex]; 1728 for (j=start; j<end; j++) { 1729 id = rbuf1_i[j] - rstart; 1730 ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 1731 sbuf2_i[j] = ncols; 1732 req_size[idex] += ncols; 1733 } 1734 req_source[idex] = r_status1[i].MPI_SOURCE; 1735 /* form the header */ 1736 sbuf2_i[0] = req_size[idex]; 1737 for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j]; 1738 ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr); 1739 } 1740 } 1741 ierr = PetscFree(r_status1);CHKERRQ(ierr); 1742 ierr = PetscFree(r_waits1);CHKERRQ(ierr); 1743 1744 /* recv buffer sizes */ 1745 /* Receive messages*/ 1746 ierr = PetscMalloc1(nrqs+1,&rbuf3);CHKERRQ(ierr); 1747 ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr); 1748 ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr); 1749 if (!ijonly) { 1750 ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr); 1751 ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr); 1752 } 1753 1754 for (i=0; i<nrqs; ++i) { 1755 ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr); 1756 ierr = PetscMalloc1(rbuf2[idex][0],&rbuf3[idex]);CHKERRQ(ierr); 1757 ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr); 1758 if (!ijonly) { 1759 ierr = PetscMalloc1(rbuf2[idex][0]*bs2,&rbuf4[idex]);CHKERRQ(ierr); 1760 ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr); 1761 } 1762 } 1763 ierr = PetscFree(r_status2);CHKERRQ(ierr); 1764 ierr = PetscFree(r_waits2);CHKERRQ(ierr); 1765 1766 /* Wait on sends1 and sends2 */ 1767 ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr); 1768 ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr); 1769 1770 if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);} 1771 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);} 1772 ierr = PetscFree(s_status1);CHKERRQ(ierr); 1773 ierr = PetscFree(s_status2);CHKERRQ(ierr); 1774 ierr = PetscFree(s_waits1);CHKERRQ(ierr); 1775 ierr = PetscFree(s_waits2);CHKERRQ(ierr); 1776 1777 /* Now allocate buffers for a->j, and send them off */ 1778 ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr); 1779 for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1780 ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr); 1781 for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 1782 1783 ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr); 1784 { 1785 for (i=0; i<nrqr; i++) { 1786 rbuf1_i = rbuf1[i]; 1787 sbuf_aj_i = sbuf_aj[i]; 1788 ct1 = 2*rbuf1_i[0] + 1; 1789 ct2 = 0; 1790 for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1791 kmax = rbuf1[i][2*j]; 1792 for (k=0; k<kmax; k++,ct1++) { 1793 row = rbuf1_i[ct1] - rstart; 1794 nzA = a_i[row+1] - a_i[row]; 1795 nzB = b_i[row+1] - b_i[row]; 1796 ncols = nzA + nzB; 1797 cworkA = a_j + a_i[row]; 1798 cworkB = b_j + b_i[row]; 1799 1800 /* load the column indices for this row into cols*/ 1801 cols = sbuf_aj_i + ct2; 1802 for (l=0; l<nzB; l++) { 1803 if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 1804 else break; 1805 } 1806 imark = l; 1807 for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 1808 for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 1809 ct2 += ncols; 1810 } 1811 } 1812 ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr); 1813 } 1814 } 1815 ierr = PetscMalloc1(nrqs+1,&r_status3);CHKERRQ(ierr); 1816 ierr = PetscMalloc1(nrqr+1,&s_status3);CHKERRQ(ierr); 1817 1818 /* Allocate buffers for a->a, and send them off */ 1819 if (!ijonly) { 1820 ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr); 1821 for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1822 ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr); 1823 for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 1824 1825 ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr); 1826 { 1827 for (i=0; i<nrqr; i++) { 1828 rbuf1_i = rbuf1[i]; 1829 sbuf_aa_i = sbuf_aa[i]; 1830 ct1 = 2*rbuf1_i[0]+1; 1831 ct2 = 0; 1832 for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1833 kmax = rbuf1_i[2*j]; 1834 for (k=0; k<kmax; k++,ct1++) { 1835 row = rbuf1_i[ct1] - rstart; 1836 nzA = a_i[row+1] - a_i[row]; 1837 nzB = b_i[row+1] - b_i[row]; 1838 ncols = nzA + nzB; 1839 cworkB = b_j + b_i[row]; 1840 vworkA = a_a + a_i[row]*bs2; 1841 vworkB = b_a + b_i[row]*bs2; 1842 1843 /* load the column values for this row into vals*/ 1844 vals = sbuf_aa_i+ct2*bs2; 1845 for (l=0; l<nzB; l++) { 1846 if ((bmap[cworkB[l]]) < cstart) { 1847 ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1848 } else break; 1849 } 1850 imark = l; 1851 for (l=0; l<nzA; l++) { 1852 ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1853 } 1854 for (l=imark; l<nzB; l++) { 1855 ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1856 } 1857 ct2 += ncols; 1858 } 1859 } 1860 ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr); 1861 } 1862 } 1863 ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr); 1864 ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr); 1865 } 1866 ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr); 1867 ierr = PetscFree(rbuf1);CHKERRQ(ierr); 1868 1869 /* Form the matrix */ 1870 /* create col map: global col of C -> local col of submatrices */ 1871 { 1872 const PetscInt *icol_i; 1873 #if defined(PETSC_USE_CTABLE) 1874 ierr = PetscMalloc1(1+ismax,&cmap);CHKERRQ(ierr); 1875 for (i=0; i<ismax; i++) { 1876 if (!allcolumns[i]) { 1877 ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr); 1878 jmax = ncol[i]; 1879 icol_i = icol[i]; 1880 cmap_i = cmap[i]; 1881 for (j=0; j<jmax; j++) { 1882 ierr = PetscTableAdd(cmap_i,icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1883 } 1884 } else { 1885 cmap[i] = NULL; 1886 } 1887 } 1888 #else 1889 ierr = PetscMalloc1(ismax,&cmap);CHKERRQ(ierr); 1890 for (i=0; i<ismax; i++) { 1891 if (!allcolumns[i]) { 1892 ierr = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr); 1893 jmax = ncol[i]; 1894 icol_i = icol[i]; 1895 cmap_i = cmap[i]; 1896 for (j=0; j<jmax; j++) { 1897 cmap_i[icol_i[j]] = j+1; 1898 } 1899 } else { /* allcolumns[i] */ 1900 cmap[i] = NULL; 1901 } 1902 } 1903 #endif 1904 } 1905 1906 /* Create lens which is required for MatCreate... */ 1907 for (i=0,j=0; i<ismax; i++) j += nrow[i]; 1908 ierr = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt),&lens);CHKERRQ(ierr); 1909 lens[0] = (PetscInt*)(lens + ismax); 1910 ierr = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr); 1911 for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1]; 1912 1913 /* Update lens from local data */ 1914 for (i=0; i<ismax; i++) { 1915 jmax = nrow[i]; 1916 if (!allcolumns[i]) cmap_i = cmap[i]; 1917 irow_i = irow[i]; 1918 lens_i = lens[i]; 1919 for (j=0; j<jmax; j++) { 1920 if (allrows[i]) row = j; 1921 else row = irow_i[j]; 1922 1923 #if defined(PETSC_USE_CTABLE) 1924 ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1925 #else 1926 proc = rtable[row]; 1927 #endif 1928 if (proc == rank) { 1929 /* Get indices from matA and then from matB */ 1930 row = row - rstart; 1931 nzA = a_i[row+1] - a_i[row]; 1932 nzB = b_i[row+1] - b_i[row]; 1933 cworkA = a_j + a_i[row]; 1934 cworkB = b_j + b_i[row]; 1935 if (!allcolumns[i]) { 1936 #if defined(PETSC_USE_CTABLE) 1937 for (k=0; k<nzA; k++) { 1938 ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1939 if (tt) lens_i[j]++; 1940 } 1941 for (k=0; k<nzB; k++) { 1942 ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1943 if (tt) lens_i[j]++; 1944 } 1945 1946 #else 1947 for (k=0; k<nzA; k++) { 1948 if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1949 } 1950 for (k=0; k<nzB; k++) { 1951 if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1952 } 1953 #endif 1954 } else { /* allcolumns */ 1955 lens_i[j] = nzA + nzB; 1956 } 1957 } 1958 } 1959 } 1960 #if defined(PETSC_USE_CTABLE) 1961 /* Create row map*/ 1962 ierr = PetscMalloc1(1+ismax,&rmap);CHKERRQ(ierr); 1963 for (i=0; i<ismax; i++) { 1964 ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr); 1965 } 1966 #else 1967 /* Create row map*/ 1968 ierr = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt),&rmap);CHKERRQ(ierr); 1969 rmap[0] = (PetscInt*)(rmap + ismax); 1970 ierr = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr); 1971 for (i=1; i<ismax; i++) rmap[i] = rmap[i-1] + Mbs; 1972 #endif 1973 for (i=0; i<ismax; i++) { 1974 irow_i = irow[i]; 1975 jmax = nrow[i]; 1976 #if defined(PETSC_USE_CTABLE) 1977 rmap_i = rmap[i]; 1978 for (j=0; j<jmax; j++) { 1979 if (allrows[i]) { 1980 ierr = PetscTableAdd(rmap_i,j+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1981 } else { 1982 ierr = PetscTableAdd(rmap_i,irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1983 } 1984 } 1985 #else 1986 rmap_i = rmap[i]; 1987 for (j=0; j<jmax; j++) { 1988 if (allrows[i]) rmap_i[j] = j; 1989 else rmap_i[irow_i[j]] = j; 1990 } 1991 #endif 1992 } 1993 1994 /* Update lens from offproc data */ 1995 { 1996 PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 1997 PetscMPIInt ii; 1998 1999 for (tmp2=0; tmp2<nrqs; tmp2++) { 2000 ierr = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr); 2001 idex = pa[ii]; 2002 sbuf1_i = sbuf1[idex]; 2003 jmax = sbuf1_i[0]; 2004 ct1 = 2*jmax+1; 2005 ct2 = 0; 2006 rbuf2_i = rbuf2[ii]; 2007 rbuf3_i = rbuf3[ii]; 2008 for (j=1; j<=jmax; j++) { 2009 is_no = sbuf1_i[2*j-1]; 2010 max1 = sbuf1_i[2*j]; 2011 lens_i = lens[is_no]; 2012 if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 2013 rmap_i = rmap[is_no]; 2014 for (k=0; k<max1; k++,ct1++) { 2015 #if defined(PETSC_USE_CTABLE) 2016 ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 2017 row--; 2018 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 2019 #else 2020 row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 2021 #endif 2022 max2 = rbuf2_i[ct1]; 2023 for (l=0; l<max2; l++,ct2++) { 2024 if (!allcolumns[is_no]) { 2025 #if defined(PETSC_USE_CTABLE) 2026 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr); 2027 if (tt) lens_i[row]++; 2028 #else 2029 if (cmap_i[rbuf3_i[ct2]]) lens_i[row]++; 2030 #endif 2031 } else { /* allcolumns */ 2032 lens_i[row]++; 2033 } 2034 } 2035 } 2036 } 2037 } 2038 } 2039 ierr = PetscFree(r_status3);CHKERRQ(ierr); 2040 ierr = PetscFree(r_waits3);CHKERRQ(ierr); 2041 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);} 2042 ierr = PetscFree(s_status3);CHKERRQ(ierr); 2043 ierr = PetscFree(s_waits3);CHKERRQ(ierr); 2044 2045 /* Create the submatrices */ 2046 if (scall == MAT_REUSE_MATRIX) { 2047 if (ijonly) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_REUSE_MATRIX and ijonly is not supported yet"); 2048 /* 2049 Assumes new rows are same length as the old rows, hence bug! 2050 */ 2051 for (i=0; i<ismax; i++) { 2052 mat = (Mat_SeqBAIJ*)(submats[i]->data); 2053 if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || C->rmap->bs != bs)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 2054 ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr); 2055 if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros"); 2056 /* Initial matrix as if empty */ 2057 ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr); 2058 2059 submats[i]->factortype = C->factortype; 2060 } 2061 } else { 2062 PetscInt bs_tmp; 2063 if (ijonly) bs_tmp = 1; 2064 else bs_tmp = bs; 2065 for (i=0; i<ismax; i++) { 2066 ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 2067 ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,nrow[i]*bs_tmp,ncol[i]*bs_tmp);CHKERRQ(ierr); 2068 ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 2069 ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); 2070 ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */ 2071 } 2072 } 2073 2074 /* Assemble the matrices */ 2075 /* First assemble the local rows */ 2076 { 2077 PetscInt ilen_row,*imat_ilen,*imat_j,*imat_i; 2078 MatScalar *imat_a = NULL; 2079 2080 for (i=0; i<ismax; i++) { 2081 mat = (Mat_SeqBAIJ*)submats[i]->data; 2082 imat_ilen = mat->ilen; 2083 imat_j = mat->j; 2084 imat_i = mat->i; 2085 if (!ijonly) imat_a = mat->a; 2086 if (!allcolumns[i]) cmap_i = cmap[i]; 2087 rmap_i = rmap[i]; 2088 irow_i = irow[i]; 2089 jmax = nrow[i]; 2090 for (j=0; j<jmax; j++) { 2091 if (allrows[i]) row = j; 2092 else row = irow_i[j]; 2093 2094 #if defined(PETSC_USE_CTABLE) 2095 ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 2096 #else 2097 proc = rtable[row]; 2098 #endif 2099 if (proc == rank) { 2100 row = row - rstart; 2101 nzA = a_i[row+1] - a_i[row]; 2102 nzB = b_i[row+1] - b_i[row]; 2103 cworkA = a_j + a_i[row]; 2104 cworkB = b_j + b_i[row]; 2105 if (!ijonly) { 2106 vworkA = a_a + a_i[row]*bs2; 2107 vworkB = b_a + b_i[row]*bs2; 2108 } 2109 #if defined(PETSC_USE_CTABLE) 2110 ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr); 2111 row--; 2112 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 2113 #else 2114 row = rmap_i[row + rstart]; 2115 #endif 2116 mat_i = imat_i[row]; 2117 if (!ijonly) mat_a = imat_a + mat_i*bs2; 2118 mat_j = imat_j + mat_i; 2119 ilen_row = imat_ilen[row]; 2120 2121 /* load the column indices for this row into cols*/ 2122 if (!allcolumns[i]) { 2123 for (l=0; l<nzB; l++) { 2124 if ((ctmp = bmap[cworkB[l]]) < cstart) { 2125 #if defined(PETSC_USE_CTABLE) 2126 ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr); 2127 if (tcol) { 2128 #else 2129 if ((tcol = cmap_i[ctmp])) { 2130 #endif 2131 *mat_j++ = tcol - 1; 2132 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2133 mat_a += bs2; 2134 ilen_row++; 2135 } 2136 } else break; 2137 } 2138 imark = l; 2139 for (l=0; l<nzA; l++) { 2140 #if defined(PETSC_USE_CTABLE) 2141 ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 2142 if (tcol) { 2143 #else 2144 if ((tcol = cmap_i[cstart + cworkA[l]])) { 2145 #endif 2146 *mat_j++ = tcol - 1; 2147 if (!ijonly) { 2148 ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2149 mat_a += bs2; 2150 } 2151 ilen_row++; 2152 } 2153 } 2154 for (l=imark; l<nzB; l++) { 2155 #if defined(PETSC_USE_CTABLE) 2156 ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 2157 if (tcol) { 2158 #else 2159 if ((tcol = cmap_i[bmap[cworkB[l]]])) { 2160 #endif 2161 *mat_j++ = tcol - 1; 2162 if (!ijonly) { 2163 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2164 mat_a += bs2; 2165 } 2166 ilen_row++; 2167 } 2168 } 2169 } else { /* allcolumns */ 2170 for (l=0; l<nzB; l++) { 2171 if ((ctmp = bmap[cworkB[l]]) < cstart) { 2172 *mat_j++ = ctmp; 2173 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2174 mat_a += bs2; 2175 ilen_row++; 2176 } else break; 2177 } 2178 imark = l; 2179 for (l=0; l<nzA; l++) { 2180 *mat_j++ = cstart+cworkA[l]; 2181 if (!ijonly) { 2182 ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2183 mat_a += bs2; 2184 } 2185 ilen_row++; 2186 } 2187 for (l=imark; l<nzB; l++) { 2188 *mat_j++ = bmap[cworkB[l]]; 2189 if (!ijonly) { 2190 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2191 mat_a += bs2; 2192 } 2193 ilen_row++; 2194 } 2195 } 2196 imat_ilen[row] = ilen_row; 2197 } 2198 } 2199 } 2200 } 2201 2202 /* Now assemble the off proc rows*/ 2203 { 2204 PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 2205 PetscInt *imat_j,*imat_i; 2206 MatScalar *imat_a = NULL,*rbuf4_i = NULL; 2207 PetscMPIInt ii; 2208 2209 for (tmp2=0; tmp2<nrqs; tmp2++) { 2210 if (ijonly) ii = tmp2; 2211 else { 2212 ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr); 2213 } 2214 idex = pa[ii]; 2215 sbuf1_i = sbuf1[idex]; 2216 jmax = sbuf1_i[0]; 2217 ct1 = 2*jmax + 1; 2218 ct2 = 0; 2219 rbuf2_i = rbuf2[ii]; 2220 rbuf3_i = rbuf3[ii]; 2221 if (!ijonly) rbuf4_i = rbuf4[ii]; 2222 for (j=1; j<=jmax; j++) { 2223 is_no = sbuf1_i[2*j-1]; 2224 if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 2225 rmap_i = rmap[is_no]; 2226 mat = (Mat_SeqBAIJ*)submats[is_no]->data; 2227 imat_ilen = mat->ilen; 2228 imat_j = mat->j; 2229 imat_i = mat->i; 2230 if (!ijonly) imat_a = mat->a; 2231 max1 = sbuf1_i[2*j]; 2232 for (k=0; k<max1; k++,ct1++) { 2233 row = sbuf1_i[ct1]; 2234 #if defined(PETSC_USE_CTABLE) 2235 ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr); 2236 row--; 2237 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 2238 #else 2239 row = rmap_i[row]; 2240 #endif 2241 ilen = imat_ilen[row]; 2242 mat_i = imat_i[row]; 2243 if (!ijonly) mat_a = imat_a + mat_i*bs2; 2244 mat_j = imat_j + mat_i; 2245 max2 = rbuf2_i[ct1]; 2246 2247 if (!allcolumns[is_no]) { 2248 for (l=0; l<max2; l++,ct2++) { 2249 #if defined(PETSC_USE_CTABLE) 2250 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 2251 if (tcol) { 2252 #else 2253 if ((tcol = cmap_i[rbuf3_i[ct2]])) { 2254 #endif 2255 *mat_j++ = tcol - 1; 2256 if (!ijonly) { 2257 ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2258 mat_a += bs2; 2259 } 2260 ilen++; 2261 } 2262 } 2263 } else { /* allcolumns */ 2264 for (l=0; l<max2; l++,ct2++) { 2265 *mat_j++ = rbuf3_i[ct2]; 2266 if (!ijonly) { 2267 ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2268 mat_a += bs2; 2269 } 2270 ilen++; 2271 } 2272 } 2273 imat_ilen[row] = ilen; 2274 } 2275 } 2276 } 2277 } 2278 /* sort the rows */ 2279 { 2280 PetscInt ilen_row,*imat_ilen,*imat_j,*imat_i; 2281 MatScalar *imat_a = NULL; 2282 MatScalar *work; 2283 2284 ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr); 2285 for (i=0; i<ismax; i++) { 2286 mat = (Mat_SeqBAIJ*)submats[i]->data; 2287 imat_ilen = mat->ilen; 2288 imat_j = mat->j; 2289 imat_i = mat->i; 2290 if (!ijonly) imat_a = mat->a; 2291 if (allcolumns[i]) continue; 2292 jmax = nrow[i]; 2293 for (j=0; j<jmax; j++) { 2294 mat_i = imat_i[j]; 2295 if (!ijonly) mat_a = imat_a + mat_i*bs2; 2296 mat_j = imat_j + mat_i; 2297 ilen_row = imat_ilen[j]; 2298 if (!ijonly) {ierr = PetscSortIntWithDataArray(ilen_row,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr);} 2299 else {ierr = PetscSortInt(ilen_row,mat_j);CHKERRQ(ierr);} 2300 } 2301 } 2302 ierr = PetscFree(work);CHKERRQ(ierr); 2303 } 2304 if (!ijonly) { 2305 ierr = PetscFree(r_status4);CHKERRQ(ierr); 2306 ierr = PetscFree(r_waits4);CHKERRQ(ierr); 2307 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);} 2308 ierr = PetscFree(s_waits4);CHKERRQ(ierr); 2309 ierr = PetscFree(s_status4);CHKERRQ(ierr); 2310 } 2311 2312 /* Restore the indices */ 2313 for (i=0; i<ismax; i++) { 2314 if (!allrows[i]) { 2315 ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 2316 } 2317 if (!allcolumns[i]) { 2318 ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 2319 } 2320 } 2321 2322 /* Destroy allocated memory */ 2323 #if defined(PETSC_USE_CTABLE) 2324 ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr); 2325 #else 2326 ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr); 2327 #endif 2328 ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 2329 ierr = PetscFree(pa);CHKERRQ(ierr); 2330 2331 ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr); 2332 ierr = PetscFree(sbuf1);CHKERRQ(ierr); 2333 ierr = PetscFree(rbuf2);CHKERRQ(ierr); 2334 for (i=0; i<nrqr; ++i) { 2335 ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 2336 } 2337 for (i=0; i<nrqs; ++i) { 2338 ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 2339 } 2340 ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr); 2341 ierr = PetscFree(rbuf3);CHKERRQ(ierr); 2342 ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 2343 ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 2344 if (!ijonly) { 2345 for (i=0; i<nrqs; ++i) {ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);} 2346 ierr = PetscFree(rbuf4);CHKERRQ(ierr); 2347 ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 2348 ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 2349 } 2350 2351 #if defined(PETSC_USE_CTABLE) 2352 for (i=0; i<ismax; i++) { 2353 ierr = PetscTableDestroy((PetscTable*)&rmap[i]);CHKERRQ(ierr); 2354 } 2355 #endif 2356 ierr = PetscFree(rmap);CHKERRQ(ierr); 2357 2358 for (i=0; i<ismax; i++) { 2359 if (!allcolumns[i]) { 2360 #if defined(PETSC_USE_CTABLE) 2361 ierr = PetscTableDestroy((PetscTable*)&cmap[i]);CHKERRQ(ierr); 2362 #else 2363 ierr = PetscFree(cmap[i]);CHKERRQ(ierr); 2364 #endif 2365 } 2366 } 2367 ierr = PetscFree(cmap);CHKERRQ(ierr); 2368 ierr = PetscFree(lens);CHKERRQ(ierr); 2369 2370 for (i=0; i<ismax; i++) { 2371 ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2372 ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2373 } 2374 2375 c->ijonly = PETSC_FALSE; /* set back to the default */ 2376 PetscFunctionReturn(0); 2377 } 2378 2379