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