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