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