1 /*$Id: mpidense.c,v 1.159 2001/08/10 03:30:41 bsmith Exp $*/ 2 3 /* 4 Basic functions for basic parallel dense matrices. 5 */ 6 7 #include "src/mat/impls/dense/mpi/mpidense.h" 8 9 EXTERN_C_BEGIN 10 #undef __FUNCT__ 11 #define __FUNCT__ "MatGetDiagonalBlock_MPIDense" 12 int MatGetDiagonalBlock_MPIDense(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *B) 13 { 14 Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 15 int m = A->m,rstart = mdn->rstart,ierr; 16 PetscScalar *array; 17 MPI_Comm comm; 18 19 PetscFunctionBegin; 20 if (A->M != A->N) SETERRQ(PETSC_ERR_SUP,"Only square matrices supported."); 21 22 /* The reuse aspect is not implemented efficiently */ 23 if (reuse) { ierr = MatDestroy(*B);CHKERRQ(ierr);} 24 25 ierr = PetscObjectGetComm((PetscObject)(mdn->A),&comm);CHKERRQ(ierr); 26 ierr = MatGetArray(mdn->A,&array);CHKERRQ(ierr); 27 ierr = MatCreate(comm,m,m,m,m,B);CHKERRQ(ierr); 28 ierr = MatSetType(*B,mdn->A->type_name);CHKERRQ(ierr); 29 ierr = MatSeqDenseSetPreallocation(*B,array+m*rstart);CHKERRQ(ierr); 30 ierr = MatRestoreArray(mdn->A,&array);CHKERRQ(ierr); 31 ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 32 ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 33 34 *iscopy = PETSC_TRUE; 35 PetscFunctionReturn(0); 36 } 37 EXTERN_C_END 38 39 #undef __FUNCT__ 40 #define __FUNCT__ "MatSetValues_MPIDense" 41 int MatSetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],const PetscScalar v[],InsertMode addv) 42 { 43 Mat_MPIDense *A = (Mat_MPIDense*)mat->data; 44 int ierr,i,j,rstart = A->rstart,rend = A->rend,row; 45 PetscTruth roworiented = A->roworiented; 46 47 PetscFunctionBegin; 48 for (i=0; i<m; i++) { 49 if (idxm[i] < 0) continue; 50 if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 51 if (idxm[i] >= rstart && idxm[i] < rend) { 52 row = idxm[i] - rstart; 53 if (roworiented) { 54 ierr = MatSetValues(A->A,1,&row,n,idxn,v+i*n,addv);CHKERRQ(ierr); 55 } else { 56 for (j=0; j<n; j++) { 57 if (idxn[j] < 0) continue; 58 if (idxn[j] >= mat->N) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 59 ierr = MatSetValues(A->A,1,&row,1,&idxn[j],v+i+j*m,addv);CHKERRQ(ierr); 60 } 61 } 62 } else { 63 if (!A->donotstash) { 64 if (roworiented) { 65 ierr = MatStashValuesRow_Private(&mat->stash,idxm[i],n,idxn,v+i*n);CHKERRQ(ierr); 66 } else { 67 ierr = MatStashValuesCol_Private(&mat->stash,idxm[i],n,idxn,v+i,m);CHKERRQ(ierr); 68 } 69 } 70 } 71 } 72 PetscFunctionReturn(0); 73 } 74 75 #undef __FUNCT__ 76 #define __FUNCT__ "MatGetValues_MPIDense" 77 int MatGetValues_MPIDense(Mat mat,int m,const int idxm[],int n,const int idxn[],PetscScalar v[]) 78 { 79 Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 80 int ierr,i,j,rstart = mdn->rstart,rend = mdn->rend,row; 81 82 PetscFunctionBegin; 83 for (i=0; i<m; i++) { 84 if (idxm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 85 if (idxm[i] >= mat->M) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 86 if (idxm[i] >= rstart && idxm[i] < rend) { 87 row = idxm[i] - rstart; 88 for (j=0; j<n; j++) { 89 if (idxn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 90 if (idxn[j] >= mat->N) { 91 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 92 } 93 ierr = MatGetValues(mdn->A,1,&row,1,&idxn[j],v+i*n+j);CHKERRQ(ierr); 94 } 95 } else { 96 SETERRQ(PETSC_ERR_SUP,"Only local values currently supported"); 97 } 98 } 99 PetscFunctionReturn(0); 100 } 101 102 #undef __FUNCT__ 103 #define __FUNCT__ "MatGetArray_MPIDense" 104 int MatGetArray_MPIDense(Mat A,PetscScalar *array[]) 105 { 106 Mat_MPIDense *a = (Mat_MPIDense*)A->data; 107 int ierr; 108 109 PetscFunctionBegin; 110 ierr = MatGetArray(a->A,array);CHKERRQ(ierr); 111 PetscFunctionReturn(0); 112 } 113 114 #undef __FUNCT__ 115 #define __FUNCT__ "MatGetSubMatrix_MPIDense" 116 static int MatGetSubMatrix_MPIDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 117 { 118 Mat_MPIDense *mat = (Mat_MPIDense*)A->data,*newmatd; 119 Mat_SeqDense *lmat = (Mat_SeqDense*)mat->A->data; 120 int i,j,ierr,*irow,*icol,rstart,rend,nrows,ncols,nlrows,nlcols; 121 PetscScalar *av,*bv,*v = lmat->v; 122 Mat newmat; 123 124 PetscFunctionBegin; 125 ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 126 ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 127 ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 128 ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 129 130 /* No parallel redistribution currently supported! Should really check each index set 131 to comfirm that it is OK. ... Currently supports only submatrix same partitioning as 132 original matrix! */ 133 134 ierr = MatGetLocalSize(A,&nlrows,&nlcols);CHKERRQ(ierr); 135 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 136 137 /* Check submatrix call */ 138 if (scall == MAT_REUSE_MATRIX) { 139 /* SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); */ 140 /* Really need to test rows and column sizes! */ 141 newmat = *B; 142 } else { 143 /* Create and fill new matrix */ 144 ierr = MatCreate(A->comm,nrows,cs,PETSC_DECIDE,ncols,&newmat);CHKERRQ(ierr); 145 ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 146 ierr = MatMPIDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 147 } 148 149 /* Now extract the data pointers and do the copy, column at a time */ 150 newmatd = (Mat_MPIDense*)newmat->data; 151 bv = ((Mat_SeqDense *)newmatd->A->data)->v; 152 153 for (i=0; i<ncols; i++) { 154 av = v + nlrows*icol[i]; 155 for (j=0; j<nrows; j++) { 156 *bv++ = av[irow[j] - rstart]; 157 } 158 } 159 160 /* Assemble the matrices so that the correct flags are set */ 161 ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 162 ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 163 164 /* Free work space */ 165 ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 166 ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 167 *B = newmat; 168 PetscFunctionReturn(0); 169 } 170 171 #undef __FUNCT__ 172 #define __FUNCT__ "MatRestoreArray_MPIDense" 173 int MatRestoreArray_MPIDense(Mat A,PetscScalar *array[]) 174 { 175 PetscFunctionBegin; 176 PetscFunctionReturn(0); 177 } 178 179 #undef __FUNCT__ 180 #define __FUNCT__ "MatAssemblyBegin_MPIDense" 181 int MatAssemblyBegin_MPIDense(Mat mat,MatAssemblyType mode) 182 { 183 Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 184 MPI_Comm comm = mat->comm; 185 int ierr,nstash,reallocs; 186 InsertMode addv; 187 188 PetscFunctionBegin; 189 /* make sure all processors are either in INSERTMODE or ADDMODE */ 190 ierr = MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);CHKERRQ(ierr); 191 if (addv == (ADD_VALUES|INSERT_VALUES)) { 192 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix adds/inserts on different procs"); 193 } 194 mat->insertmode = addv; /* in case this processor had no cache */ 195 196 ierr = MatStashScatterBegin_Private(&mat->stash,mdn->rowners);CHKERRQ(ierr); 197 ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 198 PetscLogInfo(mdn->A,"MatAssemblyBegin_MPIDense:Stash has %d entries, uses %d mallocs.\n",nstash,reallocs); 199 PetscFunctionReturn(0); 200 } 201 202 #undef __FUNCT__ 203 #define __FUNCT__ "MatAssemblyEnd_MPIDense" 204 int MatAssemblyEnd_MPIDense(Mat mat,MatAssemblyType mode) 205 { 206 Mat_MPIDense *mdn=(Mat_MPIDense*)mat->data; 207 int i,n,ierr,*row,*col,flg,j,rstart,ncols; 208 PetscScalar *val; 209 InsertMode addv=mat->insertmode; 210 211 PetscFunctionBegin; 212 /* wait on receives */ 213 while (1) { 214 ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 215 if (!flg) break; 216 217 for (i=0; i<n;) { 218 /* Now identify the consecutive vals belonging to the same row */ 219 for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; } 220 if (j < n) ncols = j-i; 221 else ncols = n-i; 222 /* Now assemble all these values with a single function call */ 223 ierr = MatSetValues_MPIDense(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr); 224 i = j; 225 } 226 } 227 ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 228 229 ierr = MatAssemblyBegin(mdn->A,mode);CHKERRQ(ierr); 230 ierr = MatAssemblyEnd(mdn->A,mode);CHKERRQ(ierr); 231 232 if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 233 ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 234 } 235 PetscFunctionReturn(0); 236 } 237 238 #undef __FUNCT__ 239 #define __FUNCT__ "MatZeroEntries_MPIDense" 240 int MatZeroEntries_MPIDense(Mat A) 241 { 242 int ierr; 243 Mat_MPIDense *l = (Mat_MPIDense*)A->data; 244 245 PetscFunctionBegin; 246 ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 247 PetscFunctionReturn(0); 248 } 249 250 #undef __FUNCT__ 251 #define __FUNCT__ "MatGetBlockSize_MPIDense" 252 int MatGetBlockSize_MPIDense(Mat A,int *bs) 253 { 254 PetscFunctionBegin; 255 *bs = 1; 256 PetscFunctionReturn(0); 257 } 258 259 /* the code does not do the diagonal entries correctly unless the 260 matrix is square and the column and row owerships are identical. 261 This is a BUG. The only way to fix it seems to be to access 262 mdn->A and mdn->B directly and not through the MatZeroRows() 263 routine. 264 */ 265 #undef __FUNCT__ 266 #define __FUNCT__ "MatZeroRows_MPIDense" 267 int MatZeroRows_MPIDense(Mat A,IS is,const PetscScalar *diag) 268 { 269 Mat_MPIDense *l = (Mat_MPIDense*)A->data; 270 int i,ierr,N,*rows,*owners = l->rowners,size = l->size; 271 int *nprocs,j,idx,nsends; 272 int nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank; 273 int *rvalues,tag = A->tag,count,base,slen,n,*source; 274 int *lens,imdex,*lrows,*values; 275 MPI_Comm comm = A->comm; 276 MPI_Request *send_waits,*recv_waits; 277 MPI_Status recv_status,*send_status; 278 IS istmp; 279 PetscTruth found; 280 281 PetscFunctionBegin; 282 ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 283 ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 284 285 /* first count number of contributors to each processor */ 286 ierr = PetscMalloc(2*size*sizeof(int),&nprocs);CHKERRQ(ierr); 287 ierr = PetscMemzero(nprocs,2*size*sizeof(int));CHKERRQ(ierr); 288 ierr = PetscMalloc((N+1)*sizeof(int),&owner);CHKERRQ(ierr); /* see note*/ 289 for (i=0; i<N; i++) { 290 idx = rows[i]; 291 found = PETSC_FALSE; 292 for (j=0; j<size; j++) { 293 if (idx >= owners[j] && idx < owners[j+1]) { 294 nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break; 295 } 296 } 297 if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index out of range"); 298 } 299 nsends = 0; for (i=0; i<size; i++) { nsends += nprocs[2*i+1];} 300 301 /* inform other processors of number of messages and max length*/ 302 ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 303 304 /* post receives: */ 305 ierr = PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);CHKERRQ(ierr); 306 ierr = PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);CHKERRQ(ierr); 307 for (i=0; i<nrecvs; i++) { 308 ierr = MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);CHKERRQ(ierr); 309 } 310 311 /* do sends: 312 1) starts[i] gives the starting index in svalues for stuff going to 313 the ith processor 314 */ 315 ierr = PetscMalloc((N+1)*sizeof(int),&svalues);CHKERRQ(ierr); 316 ierr = PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);CHKERRQ(ierr); 317 ierr = PetscMalloc((size+1)*sizeof(int),&starts);CHKERRQ(ierr); 318 starts[0] = 0; 319 for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 320 for (i=0; i<N; i++) { 321 svalues[starts[owner[i]]++] = rows[i]; 322 } 323 ISRestoreIndices(is,&rows); 324 325 starts[0] = 0; 326 for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];} 327 count = 0; 328 for (i=0; i<size; i++) { 329 if (nprocs[2*i+1]) { 330 ierr = MPI_Isend(svalues+starts[i],nprocs[2*i],MPI_INT,i,tag,comm,send_waits+count++);CHKERRQ(ierr); 331 } 332 } 333 ierr = PetscFree(starts);CHKERRQ(ierr); 334 335 base = owners[rank]; 336 337 /* wait on receives */ 338 ierr = PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);CHKERRQ(ierr); 339 source = lens + nrecvs; 340 count = nrecvs; slen = 0; 341 while (count) { 342 ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 343 /* unpack receives into our local space */ 344 ierr = MPI_Get_count(&recv_status,MPI_INT,&n);CHKERRQ(ierr); 345 source[imdex] = recv_status.MPI_SOURCE; 346 lens[imdex] = n; 347 slen += n; 348 count--; 349 } 350 ierr = PetscFree(recv_waits);CHKERRQ(ierr); 351 352 /* move the data into the send scatter */ 353 ierr = PetscMalloc((slen+1)*sizeof(int),&lrows);CHKERRQ(ierr); 354 count = 0; 355 for (i=0; i<nrecvs; i++) { 356 values = rvalues + i*nmax; 357 for (j=0; j<lens[i]; j++) { 358 lrows[count++] = values[j] - base; 359 } 360 } 361 ierr = PetscFree(rvalues);CHKERRQ(ierr); 362 ierr = PetscFree(lens);CHKERRQ(ierr); 363 ierr = PetscFree(owner);CHKERRQ(ierr); 364 ierr = PetscFree(nprocs);CHKERRQ(ierr); 365 366 /* actually zap the local rows */ 367 ierr = ISCreateGeneral(PETSC_COMM_SELF,slen,lrows,&istmp);CHKERRQ(ierr); 368 PetscLogObjectParent(A,istmp); 369 ierr = PetscFree(lrows);CHKERRQ(ierr); 370 ierr = MatZeroRows(l->A,istmp,diag);CHKERRQ(ierr); 371 ierr = ISDestroy(istmp);CHKERRQ(ierr); 372 373 /* wait on sends */ 374 if (nsends) { 375 ierr = PetscMalloc(nsends*sizeof(MPI_Status),&send_status);CHKERRQ(ierr); 376 ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 377 ierr = PetscFree(send_status);CHKERRQ(ierr); 378 } 379 ierr = PetscFree(send_waits);CHKERRQ(ierr); 380 ierr = PetscFree(svalues);CHKERRQ(ierr); 381 382 PetscFunctionReturn(0); 383 } 384 385 #undef __FUNCT__ 386 #define __FUNCT__ "MatMult_MPIDense" 387 int MatMult_MPIDense(Mat mat,Vec xx,Vec yy) 388 { 389 Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 390 int ierr; 391 392 PetscFunctionBegin; 393 ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 394 ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 395 ierr = MatMult_SeqDense(mdn->A,mdn->lvec,yy);CHKERRQ(ierr); 396 PetscFunctionReturn(0); 397 } 398 399 #undef __FUNCT__ 400 #define __FUNCT__ "MatMultAdd_MPIDense" 401 int MatMultAdd_MPIDense(Mat mat,Vec xx,Vec yy,Vec zz) 402 { 403 Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 404 int ierr; 405 406 PetscFunctionBegin; 407 ierr = VecScatterBegin(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 408 ierr = VecScatterEnd(xx,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 409 ierr = MatMultAdd_SeqDense(mdn->A,mdn->lvec,yy,zz);CHKERRQ(ierr); 410 PetscFunctionReturn(0); 411 } 412 413 #undef __FUNCT__ 414 #define __FUNCT__ "MatMultTranspose_MPIDense" 415 int MatMultTranspose_MPIDense(Mat A,Vec xx,Vec yy) 416 { 417 Mat_MPIDense *a = (Mat_MPIDense*)A->data; 418 int ierr; 419 PetscScalar zero = 0.0; 420 421 PetscFunctionBegin; 422 ierr = VecSet(&zero,yy);CHKERRQ(ierr); 423 ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 424 ierr = VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 425 ierr = VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 426 PetscFunctionReturn(0); 427 } 428 429 #undef __FUNCT__ 430 #define __FUNCT__ "MatMultTransposeAdd_MPIDense" 431 int MatMultTransposeAdd_MPIDense(Mat A,Vec xx,Vec yy,Vec zz) 432 { 433 Mat_MPIDense *a = (Mat_MPIDense*)A->data; 434 int ierr; 435 436 PetscFunctionBegin; 437 ierr = VecCopy(yy,zz);CHKERRQ(ierr); 438 ierr = MatMultTranspose_SeqDense(a->A,xx,a->lvec);CHKERRQ(ierr); 439 ierr = VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 440 ierr = VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);CHKERRQ(ierr); 441 PetscFunctionReturn(0); 442 } 443 444 #undef __FUNCT__ 445 #define __FUNCT__ "MatGetDiagonal_MPIDense" 446 int MatGetDiagonal_MPIDense(Mat A,Vec v) 447 { 448 Mat_MPIDense *a = (Mat_MPIDense*)A->data; 449 Mat_SeqDense *aloc = (Mat_SeqDense*)a->A->data; 450 int ierr,len,i,n,m = A->m,radd; 451 PetscScalar *x,zero = 0.0; 452 453 PetscFunctionBegin; 454 ierr = VecSet(&zero,v);CHKERRQ(ierr); 455 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 456 ierr = VecGetSize(v,&n);CHKERRQ(ierr); 457 if (n != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 458 len = PetscMin(a->A->m,a->A->n); 459 radd = a->rstart*m; 460 for (i=0; i<len; i++) { 461 x[i] = aloc->v[radd + i*m + i]; 462 } 463 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 464 PetscFunctionReturn(0); 465 } 466 467 #undef __FUNCT__ 468 #define __FUNCT__ "MatDestroy_MPIDense" 469 int MatDestroy_MPIDense(Mat mat) 470 { 471 Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 472 int ierr; 473 474 PetscFunctionBegin; 475 476 #if defined(PETSC_USE_LOG) 477 PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->M,mat->N); 478 #endif 479 ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 480 ierr = PetscFree(mdn->rowners);CHKERRQ(ierr); 481 ierr = MatDestroy(mdn->A);CHKERRQ(ierr); 482 if (mdn->lvec) VecDestroy(mdn->lvec); 483 if (mdn->Mvctx) VecScatterDestroy(mdn->Mvctx); 484 if (mdn->factor) { 485 if (mdn->factor->temp) {ierr = PetscFree(mdn->factor->temp);CHKERRQ(ierr);} 486 if (mdn->factor->tag) {ierr = PetscFree(mdn->factor->tag);CHKERRQ(ierr);} 487 if (mdn->factor->pivots) {ierr = PetscFree(mdn->factor->pivots);CHKERRQ(ierr);} 488 ierr = PetscFree(mdn->factor);CHKERRQ(ierr); 489 } 490 ierr = PetscFree(mdn);CHKERRQ(ierr); 491 PetscFunctionReturn(0); 492 } 493 494 #undef __FUNCT__ 495 #define __FUNCT__ "MatView_MPIDense_Binary" 496 static int MatView_MPIDense_Binary(Mat mat,PetscViewer viewer) 497 { 498 Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 499 int ierr; 500 501 PetscFunctionBegin; 502 if (mdn->size == 1) { 503 ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 504 } 505 else SETERRQ(PETSC_ERR_SUP,"Only uniprocessor output supported"); 506 PetscFunctionReturn(0); 507 } 508 509 #undef __FUNCT__ 510 #define __FUNCT__ "MatView_MPIDense_ASCIIorDraworSocket" 511 static int MatView_MPIDense_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 512 { 513 Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data; 514 int ierr,size = mdn->size,rank = mdn->rank; 515 PetscViewerType vtype; 516 PetscTruth isascii,isdraw; 517 PetscViewer sviewer; 518 PetscViewerFormat format; 519 520 PetscFunctionBegin; 521 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 522 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 523 if (isascii) { 524 ierr = PetscViewerGetType(viewer,&vtype);CHKERRQ(ierr); 525 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 526 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 527 MatInfo info; 528 ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 529 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] local rows %d nz %d nz alloced %d mem %d \n",rank,mat->m, 530 (int)info.nz_used,(int)info.nz_allocated,(int)info.memory);CHKERRQ(ierr); 531 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 532 ierr = VecScatterView(mdn->Mvctx,viewer);CHKERRQ(ierr); 533 PetscFunctionReturn(0); 534 } else if (format == PETSC_VIEWER_ASCII_INFO) { 535 PetscFunctionReturn(0); 536 } 537 } else if (isdraw) { 538 PetscDraw draw; 539 PetscTruth isnull; 540 541 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 542 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 543 if (isnull) PetscFunctionReturn(0); 544 } 545 546 if (size == 1) { 547 ierr = MatView(mdn->A,viewer);CHKERRQ(ierr); 548 } else { 549 /* assemble the entire matrix onto first processor. */ 550 Mat A; 551 int M = mat->M,N = mat->N,m,row,i,nz,*cols; 552 PetscScalar *vals; 553 554 if (!rank) { 555 ierr = MatCreate(mat->comm,M,N,M,N,&A);CHKERRQ(ierr); 556 } else { 557 ierr = MatCreate(mat->comm,0,0,M,N,&A);CHKERRQ(ierr); 558 } 559 /* Since this is a temporary matrix, MATMPIDENSE instead of A->type_name here is probably acceptable. */ 560 ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 561 ierr = MatMPIDenseSetPreallocation(A,PETSC_NULL); 562 PetscLogObjectParent(mat,A); 563 564 /* Copy the matrix ... This isn't the most efficient means, 565 but it's quick for now */ 566 row = mdn->rstart; m = mdn->A->m; 567 for (i=0; i<m; i++) { 568 ierr = MatGetRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 569 ierr = MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 570 ierr = MatRestoreRow(mat,row,&nz,&cols,&vals);CHKERRQ(ierr); 571 row++; 572 } 573 574 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 575 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 576 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 577 if (!rank) { 578 ierr = MatView(((Mat_MPIDense*)(A->data))->A,sviewer);CHKERRQ(ierr); 579 } 580 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 581 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 582 ierr = MatDestroy(A);CHKERRQ(ierr); 583 } 584 PetscFunctionReturn(0); 585 } 586 587 #undef __FUNCT__ 588 #define __FUNCT__ "MatView_MPIDense" 589 int MatView_MPIDense(Mat mat,PetscViewer viewer) 590 { 591 int ierr; 592 PetscTruth isascii,isbinary,isdraw,issocket; 593 594 PetscFunctionBegin; 595 596 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 597 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 598 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 599 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 600 601 if (isascii || issocket || isdraw) { 602 ierr = MatView_MPIDense_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 603 } else if (isbinary) { 604 ierr = MatView_MPIDense_Binary(mat,viewer);CHKERRQ(ierr); 605 } else { 606 SETERRQ1(1,"Viewer type %s not supported by MPI dense matrix",((PetscObject)viewer)->type_name); 607 } 608 PetscFunctionReturn(0); 609 } 610 611 #undef __FUNCT__ 612 #define __FUNCT__ "MatGetInfo_MPIDense" 613 int MatGetInfo_MPIDense(Mat A,MatInfoType flag,MatInfo *info) 614 { 615 Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 616 Mat mdn = mat->A; 617 int ierr; 618 PetscReal isend[5],irecv[5]; 619 620 PetscFunctionBegin; 621 info->rows_global = (double)A->M; 622 info->columns_global = (double)A->N; 623 info->rows_local = (double)A->m; 624 info->columns_local = (double)A->N; 625 info->block_size = 1.0; 626 ierr = MatGetInfo(mdn,MAT_LOCAL,info);CHKERRQ(ierr); 627 isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 628 isend[3] = info->memory; isend[4] = info->mallocs; 629 if (flag == MAT_LOCAL) { 630 info->nz_used = isend[0]; 631 info->nz_allocated = isend[1]; 632 info->nz_unneeded = isend[2]; 633 info->memory = isend[3]; 634 info->mallocs = isend[4]; 635 } else if (flag == MAT_GLOBAL_MAX) { 636 ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 637 info->nz_used = irecv[0]; 638 info->nz_allocated = irecv[1]; 639 info->nz_unneeded = irecv[2]; 640 info->memory = irecv[3]; 641 info->mallocs = irecv[4]; 642 } else if (flag == MAT_GLOBAL_SUM) { 643 ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 644 info->nz_used = irecv[0]; 645 info->nz_allocated = irecv[1]; 646 info->nz_unneeded = irecv[2]; 647 info->memory = irecv[3]; 648 info->mallocs = irecv[4]; 649 } 650 info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 651 info->fill_ratio_needed = 0; 652 info->factor_mallocs = 0; 653 PetscFunctionReturn(0); 654 } 655 656 #undef __FUNCT__ 657 #define __FUNCT__ "MatSetOption_MPIDense" 658 int MatSetOption_MPIDense(Mat A,MatOption op) 659 { 660 Mat_MPIDense *a = (Mat_MPIDense*)A->data; 661 int ierr; 662 663 PetscFunctionBegin; 664 switch (op) { 665 case MAT_NO_NEW_NONZERO_LOCATIONS: 666 case MAT_YES_NEW_NONZERO_LOCATIONS: 667 case MAT_NEW_NONZERO_LOCATION_ERR: 668 case MAT_NEW_NONZERO_ALLOCATION_ERR: 669 case MAT_COLUMNS_SORTED: 670 case MAT_COLUMNS_UNSORTED: 671 ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 672 break; 673 case MAT_ROW_ORIENTED: 674 a->roworiented = PETSC_TRUE; 675 ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 676 break; 677 case MAT_ROWS_SORTED: 678 case MAT_ROWS_UNSORTED: 679 case MAT_YES_NEW_DIAGONALS: 680 case MAT_USE_HASH_TABLE: 681 PetscLogInfo(A,"MatSetOption_MPIDense:Option ignored\n"); 682 break; 683 case MAT_COLUMN_ORIENTED: 684 a->roworiented = PETSC_FALSE; 685 ierr = MatSetOption(a->A,op);CHKERRQ(ierr); 686 break; 687 case MAT_IGNORE_OFF_PROC_ENTRIES: 688 a->donotstash = PETSC_TRUE; 689 break; 690 case MAT_NO_NEW_DIAGONALS: 691 SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 692 case MAT_SYMMETRIC: 693 case MAT_STRUCTURALLY_SYMMETRIC: 694 case MAT_NOT_SYMMETRIC: 695 case MAT_NOT_STRUCTURALLY_SYMMETRIC: 696 case MAT_HERMITIAN: 697 case MAT_NOT_HERMITIAN: 698 case MAT_SYMMETRY_ETERNAL: 699 case MAT_NOT_SYMMETRY_ETERNAL: 700 break; 701 default: 702 SETERRQ(PETSC_ERR_SUP,"unknown option"); 703 } 704 PetscFunctionReturn(0); 705 } 706 707 #undef __FUNCT__ 708 #define __FUNCT__ "MatGetRow_MPIDense" 709 int MatGetRow_MPIDense(Mat A,int row,int *nz,int **idx,PetscScalar **v) 710 { 711 Mat_MPIDense *mat = (Mat_MPIDense*)A->data; 712 int lrow,rstart = mat->rstart,rend = mat->rend,ierr; 713 714 PetscFunctionBegin; 715 if (row < rstart || row >= rend) SETERRQ(PETSC_ERR_SUP,"only local rows") 716 lrow = row - rstart; 717 ierr = MatGetRow(mat->A,lrow,nz,idx,v);CHKERRQ(ierr); 718 PetscFunctionReturn(0); 719 } 720 721 #undef __FUNCT__ 722 #define __FUNCT__ "MatRestoreRow_MPIDense" 723 int MatRestoreRow_MPIDense(Mat mat,int row,int *nz,int **idx,PetscScalar **v) 724 { 725 int ierr; 726 727 PetscFunctionBegin; 728 if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 729 if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 730 PetscFunctionReturn(0); 731 } 732 733 #undef __FUNCT__ 734 #define __FUNCT__ "MatDiagonalScale_MPIDense" 735 int MatDiagonalScale_MPIDense(Mat A,Vec ll,Vec rr) 736 { 737 Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 738 Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 739 PetscScalar *l,*r,x,*v; 740 int ierr,i,j,s2a,s3a,s2,s3,m=mdn->A->m,n=mdn->A->n; 741 742 PetscFunctionBegin; 743 ierr = MatGetLocalSize(A,&s2,&s3);CHKERRQ(ierr); 744 if (ll) { 745 ierr = VecGetLocalSize(ll,&s2a);CHKERRQ(ierr); 746 if (s2a != s2) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector non-conforming local size"); 747 ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 748 for (i=0; i<m; i++) { 749 x = l[i]; 750 v = mat->v + i; 751 for (j=0; j<n; j++) { (*v) *= x; v+= m;} 752 } 753 ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 754 PetscLogFlops(n*m); 755 } 756 if (rr) { 757 ierr = VecGetSize(rr,&s3a);CHKERRQ(ierr); 758 if (s3a != s3) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec non-conforming local size"); 759 ierr = VecScatterBegin(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 760 ierr = VecScatterEnd(rr,mdn->lvec,INSERT_VALUES,SCATTER_FORWARD,mdn->Mvctx);CHKERRQ(ierr); 761 ierr = VecGetArray(mdn->lvec,&r);CHKERRQ(ierr); 762 for (i=0; i<n; i++) { 763 x = r[i]; 764 v = mat->v + i*m; 765 for (j=0; j<m; j++) { (*v++) *= x;} 766 } 767 ierr = VecRestoreArray(mdn->lvec,&r);CHKERRQ(ierr); 768 PetscLogFlops(n*m); 769 } 770 PetscFunctionReturn(0); 771 } 772 773 #undef __FUNCT__ 774 #define __FUNCT__ "MatNorm_MPIDense" 775 int MatNorm_MPIDense(Mat A,NormType type,PetscReal *nrm) 776 { 777 Mat_MPIDense *mdn = (Mat_MPIDense*)A->data; 778 Mat_SeqDense *mat = (Mat_SeqDense*)mdn->A->data; 779 int ierr,i,j; 780 PetscReal sum = 0.0; 781 PetscScalar *v = mat->v; 782 783 PetscFunctionBegin; 784 if (mdn->size == 1) { 785 ierr = MatNorm(mdn->A,type,nrm);CHKERRQ(ierr); 786 } else { 787 if (type == NORM_FROBENIUS) { 788 for (i=0; i<mdn->A->n*mdn->A->m; i++) { 789 #if defined(PETSC_USE_COMPLEX) 790 sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 791 #else 792 sum += (*v)*(*v); v++; 793 #endif 794 } 795 ierr = MPI_Allreduce(&sum,nrm,1,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 796 *nrm = sqrt(*nrm); 797 PetscLogFlops(2*mdn->A->n*mdn->A->m); 798 } else if (type == NORM_1) { 799 PetscReal *tmp,*tmp2; 800 ierr = PetscMalloc(2*A->N*sizeof(PetscReal),&tmp);CHKERRQ(ierr); 801 tmp2 = tmp + A->N; 802 ierr = PetscMemzero(tmp,2*A->N*sizeof(PetscReal));CHKERRQ(ierr); 803 *nrm = 0.0; 804 v = mat->v; 805 for (j=0; j<mdn->A->n; j++) { 806 for (i=0; i<mdn->A->m; i++) { 807 tmp[j] += PetscAbsScalar(*v); v++; 808 } 809 } 810 ierr = MPI_Allreduce(tmp,tmp2,A->N,MPIU_REAL,MPI_SUM,A->comm);CHKERRQ(ierr); 811 for (j=0; j<A->N; j++) { 812 if (tmp2[j] > *nrm) *nrm = tmp2[j]; 813 } 814 ierr = PetscFree(tmp);CHKERRQ(ierr); 815 PetscLogFlops(A->n*A->m); 816 } else if (type == NORM_INFINITY) { /* max row norm */ 817 PetscReal ntemp; 818 ierr = MatNorm(mdn->A,type,&ntemp);CHKERRQ(ierr); 819 ierr = MPI_Allreduce(&ntemp,nrm,1,MPIU_REAL,MPI_MAX,A->comm);CHKERRQ(ierr); 820 } else { 821 SETERRQ(PETSC_ERR_SUP,"No support for two norm"); 822 } 823 } 824 PetscFunctionReturn(0); 825 } 826 827 #undef __FUNCT__ 828 #define __FUNCT__ "MatTranspose_MPIDense" 829 int MatTranspose_MPIDense(Mat A,Mat *matout) 830 { 831 Mat_MPIDense *a = (Mat_MPIDense*)A->data; 832 Mat_SeqDense *Aloc = (Mat_SeqDense*)a->A->data; 833 Mat B; 834 int M = A->M,N = A->N,m,n,*rwork,rstart = a->rstart; 835 int j,i,ierr; 836 PetscScalar *v; 837 838 PetscFunctionBegin; 839 if (!matout && M != N) { 840 SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place"); 841 } 842 ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,N,M,&B);CHKERRQ(ierr); 843 ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 844 ierr = MatMPIDenseSetPreallocation(B,PETSC_NULL);CHKERRQ(ierr); 845 846 m = a->A->m; n = a->A->n; v = Aloc->v; 847 ierr = PetscMalloc(n*sizeof(int),&rwork);CHKERRQ(ierr); 848 for (j=0; j<n; j++) { 849 for (i=0; i<m; i++) rwork[i] = rstart + i; 850 ierr = MatSetValues(B,1,&j,m,rwork,v,INSERT_VALUES);CHKERRQ(ierr); 851 v += m; 852 } 853 ierr = PetscFree(rwork);CHKERRQ(ierr); 854 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 855 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 856 if (matout) { 857 *matout = B; 858 } else { 859 ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 860 } 861 PetscFunctionReturn(0); 862 } 863 864 #include "petscblaslapack.h" 865 #undef __FUNCT__ 866 #define __FUNCT__ "MatScale_MPIDense" 867 int MatScale_MPIDense(const PetscScalar *alpha,Mat inA) 868 { 869 Mat_MPIDense *A = (Mat_MPIDense*)inA->data; 870 Mat_SeqDense *a = (Mat_SeqDense*)A->A->data; 871 int one = 1,nz; 872 873 PetscFunctionBegin; 874 nz = inA->m*inA->N; 875 BLscal_(&nz,(PetscScalar*)alpha,a->v,&one); 876 PetscLogFlops(nz); 877 PetscFunctionReturn(0); 878 } 879 880 static int MatDuplicate_MPIDense(Mat,MatDuplicateOption,Mat *); 881 882 #undef __FUNCT__ 883 #define __FUNCT__ "MatSetUpPreallocation_MPIDense" 884 int MatSetUpPreallocation_MPIDense(Mat A) 885 { 886 int ierr; 887 888 PetscFunctionBegin; 889 ierr = MatMPIDenseSetPreallocation(A,0);CHKERRQ(ierr); 890 PetscFunctionReturn(0); 891 } 892 893 /* -------------------------------------------------------------------*/ 894 static struct _MatOps MatOps_Values = {MatSetValues_MPIDense, 895 MatGetRow_MPIDense, 896 MatRestoreRow_MPIDense, 897 MatMult_MPIDense, 898 /* 4*/ MatMultAdd_MPIDense, 899 MatMultTranspose_MPIDense, 900 MatMultTransposeAdd_MPIDense, 901 0, 902 0, 903 0, 904 /*10*/ 0, 905 0, 906 0, 907 0, 908 MatTranspose_MPIDense, 909 /*15*/ MatGetInfo_MPIDense, 910 0, 911 MatGetDiagonal_MPIDense, 912 MatDiagonalScale_MPIDense, 913 MatNorm_MPIDense, 914 /*20*/ MatAssemblyBegin_MPIDense, 915 MatAssemblyEnd_MPIDense, 916 0, 917 MatSetOption_MPIDense, 918 MatZeroEntries_MPIDense, 919 /*25*/ MatZeroRows_MPIDense, 920 0, 921 0, 922 0, 923 0, 924 /*30*/ MatSetUpPreallocation_MPIDense, 925 0, 926 0, 927 MatGetArray_MPIDense, 928 MatRestoreArray_MPIDense, 929 /*35*/ MatDuplicate_MPIDense, 930 0, 931 0, 932 0, 933 0, 934 /*40*/ 0, 935 MatGetSubMatrices_MPIDense, 936 0, 937 MatGetValues_MPIDense, 938 0, 939 /*45*/ 0, 940 MatScale_MPIDense, 941 0, 942 0, 943 0, 944 /*50*/ MatGetBlockSize_MPIDense, 945 0, 946 0, 947 0, 948 0, 949 /*55*/ 0, 950 0, 951 0, 952 0, 953 0, 954 /*60*/ MatGetSubMatrix_MPIDense, 955 MatDestroy_MPIDense, 956 MatView_MPIDense, 957 MatGetPetscMaps_Petsc, 958 0, 959 /*65*/ 0, 960 0, 961 0, 962 0, 963 0, 964 /*70*/ 0, 965 0, 966 0, 967 0, 968 0, 969 /*75*/ 0, 970 0, 971 0, 972 0, 973 0, 974 /*80*/ 0, 975 0, 976 0, 977 0, 978 /*85*/ MatLoad_MPIDense}; 979 980 EXTERN_C_BEGIN 981 #undef __FUNCT__ 982 #define __FUNCT__ "MatMPIDenseSetPreallocation_MPIDense" 983 int MatMPIDenseSetPreallocation_MPIDense(Mat mat,PetscScalar *data) 984 { 985 Mat_MPIDense *a; 986 int ierr; 987 988 PetscFunctionBegin; 989 mat->preallocated = PETSC_TRUE; 990 /* Note: For now, when data is specified above, this assumes the user correctly 991 allocates the local dense storage space. We should add error checking. */ 992 993 a = (Mat_MPIDense*)mat->data; 994 ierr = MatCreate(PETSC_COMM_SELF,mat->m,mat->N,mat->m,mat->N,&a->A);CHKERRQ(ierr); 995 ierr = MatSetType(a->A,MATSEQDENSE);CHKERRQ(ierr); 996 ierr = MatSeqDenseSetPreallocation(a->A,data);CHKERRQ(ierr); 997 PetscLogObjectParent(mat,a->A); 998 PetscFunctionReturn(0); 999 } 1000 EXTERN_C_END 1001 1002 /*MC 1003 MATMPIDENSE - MATMPIDENSE = "mpidense" - A matrix type to be used for distributed dense matrices. 1004 1005 Options Database Keys: 1006 . -mat_type mpidense - sets the matrix type to "mpidense" during a call to MatSetFromOptions() 1007 1008 Level: beginner 1009 1010 .seealso: MatCreateMPIDense 1011 M*/ 1012 1013 EXTERN_C_BEGIN 1014 #undef __FUNCT__ 1015 #define __FUNCT__ "MatCreate_MPIDense" 1016 int MatCreate_MPIDense(Mat mat) 1017 { 1018 Mat_MPIDense *a; 1019 int ierr,i; 1020 1021 PetscFunctionBegin; 1022 ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1023 mat->data = (void*)a; 1024 ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1025 mat->factor = 0; 1026 mat->mapping = 0; 1027 1028 a->factor = 0; 1029 mat->insertmode = NOT_SET_VALUES; 1030 ierr = MPI_Comm_rank(mat->comm,&a->rank);CHKERRQ(ierr); 1031 ierr = MPI_Comm_size(mat->comm,&a->size);CHKERRQ(ierr); 1032 1033 ierr = PetscSplitOwnership(mat->comm,&mat->m,&mat->M);CHKERRQ(ierr); 1034 ierr = PetscSplitOwnership(mat->comm,&mat->n,&mat->N);CHKERRQ(ierr); 1035 a->nvec = mat->n; 1036 1037 /* the information in the maps duplicates the information computed below, eventually 1038 we should remove the duplicate information that is not contained in the maps */ 1039 ierr = PetscMapCreateMPI(mat->comm,mat->m,mat->M,&mat->rmap);CHKERRQ(ierr); 1040 ierr = PetscMapCreateMPI(mat->comm,mat->n,mat->N,&mat->cmap);CHKERRQ(ierr); 1041 1042 /* build local table of row and column ownerships */ 1043 ierr = PetscMalloc(2*(a->size+2)*sizeof(int),&a->rowners);CHKERRQ(ierr); 1044 a->cowners = a->rowners + a->size + 1; 1045 PetscLogObjectMemory(mat,2*(a->size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1046 ierr = MPI_Allgather(&mat->m,1,MPI_INT,a->rowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr); 1047 a->rowners[0] = 0; 1048 for (i=2; i<=a->size; i++) { 1049 a->rowners[i] += a->rowners[i-1]; 1050 } 1051 a->rstart = a->rowners[a->rank]; 1052 a->rend = a->rowners[a->rank+1]; 1053 ierr = MPI_Allgather(&mat->n,1,MPI_INT,a->cowners+1,1,MPI_INT,mat->comm);CHKERRQ(ierr); 1054 a->cowners[0] = 0; 1055 for (i=2; i<=a->size; i++) { 1056 a->cowners[i] += a->cowners[i-1]; 1057 } 1058 1059 /* build cache for off array entries formed */ 1060 a->donotstash = PETSC_FALSE; 1061 ierr = MatStashCreate_Private(mat->comm,1,&mat->stash);CHKERRQ(ierr); 1062 1063 /* stuff used for matrix vector multiply */ 1064 a->lvec = 0; 1065 a->Mvctx = 0; 1066 a->roworiented = PETSC_TRUE; 1067 1068 ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatGetDiagonalBlock_C", 1069 "MatGetDiagonalBlock_MPIDense", 1070 MatGetDiagonalBlock_MPIDense);CHKERRQ(ierr); 1071 ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMPIDenseSetPreallocation_C", 1072 "MatMPIDenseSetPreallocation_MPIDense", 1073 MatMPIDenseSetPreallocation_MPIDense);CHKERRQ(ierr); 1074 PetscFunctionReturn(0); 1075 } 1076 EXTERN_C_END 1077 1078 /*MC 1079 MATDENSE - MATDENSE = "dense" - A matrix type to be used for dense matrices. 1080 1081 This matrix type is identical to MATSEQDENSE when constructed with a single process communicator, 1082 and MATMPIDENSE otherwise. 1083 1084 Options Database Keys: 1085 . -mat_type dense - sets the matrix type to "dense" during a call to MatSetFromOptions() 1086 1087 Level: beginner 1088 1089 .seealso: MatCreateMPIDense,MATSEQDENSE,MATMPIDENSE 1090 M*/ 1091 1092 EXTERN_C_BEGIN 1093 #undef __FUNCT__ 1094 #define __FUNCT__ "MatCreate_Dense" 1095 int MatCreate_Dense(Mat A) { 1096 int ierr,size; 1097 1098 PetscFunctionBegin; 1099 ierr = PetscObjectChangeTypeName((PetscObject)A,MATDENSE);CHKERRQ(ierr); 1100 ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 1101 if (size == 1) { 1102 ierr = MatSetType(A,MATSEQDENSE);CHKERRQ(ierr); 1103 } else { 1104 ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr); 1105 } 1106 PetscFunctionReturn(0); 1107 } 1108 EXTERN_C_END 1109 1110 #undef __FUNCT__ 1111 #define __FUNCT__ "MatMPIDenseSetPreallocation" 1112 /*@C 1113 MatMPIDenseSetPreallocation - Sets the array used to store the matrix entries 1114 1115 Not collective 1116 1117 Input Parameters: 1118 . A - the matrix 1119 - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1120 to control all matrix memory allocation. 1121 1122 Notes: 1123 The dense format is fully compatible with standard Fortran 77 1124 storage by columns. 1125 1126 The data input variable is intended primarily for Fortran programmers 1127 who wish to allocate their own matrix memory space. Most users should 1128 set data=PETSC_NULL. 1129 1130 Level: intermediate 1131 1132 .keywords: matrix,dense, parallel 1133 1134 .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1135 @*/ 1136 int MatMPIDenseSetPreallocation(Mat mat,PetscScalar *data) 1137 { 1138 int ierr,(*f)(Mat,PetscScalar *); 1139 1140 PetscFunctionBegin; 1141 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMPIDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1142 if (f) { 1143 ierr = (*f)(mat,data);CHKERRQ(ierr); 1144 } 1145 PetscFunctionReturn(0); 1146 } 1147 1148 #undef __FUNCT__ 1149 #define __FUNCT__ "MatCreateMPIDense" 1150 /*@C 1151 MatCreateMPIDense - Creates a sparse parallel matrix in dense format. 1152 1153 Collective on MPI_Comm 1154 1155 Input Parameters: 1156 + comm - MPI communicator 1157 . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1158 . n - number of local columns (or PETSC_DECIDE to have calculated if N is given) 1159 . M - number of global rows (or PETSC_DECIDE to have calculated if m is given) 1160 . N - number of global columns (or PETSC_DECIDE to have calculated if n is given) 1161 - data - optional location of matrix data. Set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users) for PETSc 1162 to control all matrix memory allocation. 1163 1164 Output Parameter: 1165 . A - the matrix 1166 1167 Notes: 1168 The dense format is fully compatible with standard Fortran 77 1169 storage by columns. 1170 1171 The data input variable is intended primarily for Fortran programmers 1172 who wish to allocate their own matrix memory space. Most users should 1173 set data=PETSC_NULL (PETSC_NULL_SCALAR for Fortran users). 1174 1175 The user MUST specify either the local or global matrix dimensions 1176 (possibly both). 1177 1178 Level: intermediate 1179 1180 .keywords: matrix,dense, parallel 1181 1182 .seealso: MatCreate(), MatCreateSeqDense(), MatSetValues() 1183 @*/ 1184 int MatCreateMPIDense(MPI_Comm comm,int m,int n,int M,int N,PetscScalar *data,Mat *A) 1185 { 1186 int ierr,size; 1187 1188 PetscFunctionBegin; 1189 ierr = MatCreate(comm,m,n,M,N,A);CHKERRQ(ierr); 1190 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1191 if (size > 1) { 1192 ierr = MatSetType(*A,MATMPIDENSE);CHKERRQ(ierr); 1193 ierr = MatMPIDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1194 } else { 1195 ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1196 ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1197 } 1198 PetscFunctionReturn(0); 1199 } 1200 1201 #undef __FUNCT__ 1202 #define __FUNCT__ "MatDuplicate_MPIDense" 1203 static int MatDuplicate_MPIDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 1204 { 1205 Mat mat; 1206 Mat_MPIDense *a,*oldmat = (Mat_MPIDense*)A->data; 1207 int ierr; 1208 1209 PetscFunctionBegin; 1210 *newmat = 0; 1211 ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&mat);CHKERRQ(ierr); 1212 ierr = MatSetType(mat,A->type_name);CHKERRQ(ierr); 1213 ierr = PetscNew(Mat_MPIDense,&a);CHKERRQ(ierr); 1214 mat->data = (void*)a; 1215 ierr = PetscMemcpy(mat->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1216 mat->factor = A->factor; 1217 mat->assembled = PETSC_TRUE; 1218 mat->preallocated = PETSC_TRUE; 1219 1220 a->rstart = oldmat->rstart; 1221 a->rend = oldmat->rend; 1222 a->size = oldmat->size; 1223 a->rank = oldmat->rank; 1224 mat->insertmode = NOT_SET_VALUES; 1225 a->nvec = oldmat->nvec; 1226 a->donotstash = oldmat->donotstash; 1227 ierr = PetscMalloc((a->size+1)*sizeof(int),&a->rowners);CHKERRQ(ierr); 1228 PetscLogObjectMemory(mat,(a->size+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIDense)); 1229 ierr = PetscMemcpy(a->rowners,oldmat->rowners,(a->size+1)*sizeof(int));CHKERRQ(ierr); 1230 ierr = MatStashCreate_Private(A->comm,1,&mat->stash);CHKERRQ(ierr); 1231 1232 ierr = MatSetUpMultiply_MPIDense(mat);CHKERRQ(ierr); 1233 ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 1234 PetscLogObjectParent(mat,a->A); 1235 *newmat = mat; 1236 PetscFunctionReturn(0); 1237 } 1238 1239 #include "petscsys.h" 1240 1241 #undef __FUNCT__ 1242 #define __FUNCT__ "MatLoad_MPIDense_DenseInFile" 1243 int MatLoad_MPIDense_DenseInFile(MPI_Comm comm,int fd,int M,int N,const MatType type,Mat *newmat) 1244 { 1245 int *rowners,i,size,rank,m,ierr,nz,j; 1246 PetscScalar *array,*vals,*vals_ptr; 1247 MPI_Status status; 1248 1249 PetscFunctionBegin; 1250 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1251 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1252 1253 /* determine ownership of all rows */ 1254 m = M/size + ((M % size) > rank); 1255 ierr = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr); 1256 ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 1257 rowners[0] = 0; 1258 for (i=2; i<=size; i++) { 1259 rowners[i] += rowners[i-1]; 1260 } 1261 1262 ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr); 1263 ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1264 ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 1265 ierr = MatGetArray(*newmat,&array);CHKERRQ(ierr); 1266 1267 if (!rank) { 1268 ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 1269 1270 /* read in my part of the matrix numerical values */ 1271 ierr = PetscBinaryRead(fd,vals,m*N,PETSC_SCALAR);CHKERRQ(ierr); 1272 1273 /* insert into matrix-by row (this is why cannot directly read into array */ 1274 vals_ptr = vals; 1275 for (i=0; i<m; i++) { 1276 for (j=0; j<N; j++) { 1277 array[i + j*m] = *vals_ptr++; 1278 } 1279 } 1280 1281 /* read in other processors and ship out */ 1282 for (i=1; i<size; i++) { 1283 nz = (rowners[i+1] - rowners[i])*N; 1284 ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1285 ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,(*newmat)->tag,comm);CHKERRQ(ierr); 1286 } 1287 } else { 1288 /* receive numeric values */ 1289 ierr = PetscMalloc(m*N*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 1290 1291 /* receive message of values*/ 1292 ierr = MPI_Recv(vals,m*N,MPIU_SCALAR,0,(*newmat)->tag,comm,&status);CHKERRQ(ierr); 1293 1294 /* insert into matrix-by row (this is why cannot directly read into array */ 1295 vals_ptr = vals; 1296 for (i=0; i<m; i++) { 1297 for (j=0; j<N; j++) { 1298 array[i + j*m] = *vals_ptr++; 1299 } 1300 } 1301 } 1302 ierr = PetscFree(rowners);CHKERRQ(ierr); 1303 ierr = PetscFree(vals);CHKERRQ(ierr); 1304 ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1305 ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1306 PetscFunctionReturn(0); 1307 } 1308 1309 #undef __FUNCT__ 1310 #define __FUNCT__ "MatLoad_MPIDense" 1311 int MatLoad_MPIDense(PetscViewer viewer,const MatType type,Mat *newmat) 1312 { 1313 Mat A; 1314 PetscScalar *vals,*svals; 1315 MPI_Comm comm = ((PetscObject)viewer)->comm; 1316 MPI_Status status; 1317 int header[4],rank,size,*rowlengths = 0,M,N,m,*rowners,maxnz,*cols; 1318 int *ourlens,*sndcounts = 0,*procsnz = 0,*offlens,jj,*mycols,*smycols; 1319 int tag = ((PetscObject)viewer)->tag; 1320 int i,nz,ierr,j,rstart,rend,fd; 1321 1322 PetscFunctionBegin; 1323 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1324 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1325 if (!rank) { 1326 ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 1327 ierr = PetscBinaryRead(fd,(char *)header,4,PETSC_INT);CHKERRQ(ierr); 1328 if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object"); 1329 } 1330 1331 ierr = MPI_Bcast(header+1,3,MPI_INT,0,comm);CHKERRQ(ierr); 1332 M = header[1]; N = header[2]; nz = header[3]; 1333 1334 /* 1335 Handle case where matrix is stored on disk as a dense matrix 1336 */ 1337 if (nz == MATRIX_BINARY_FORMAT_DENSE) { 1338 ierr = MatLoad_MPIDense_DenseInFile(comm,fd,M,N,type,newmat);CHKERRQ(ierr); 1339 PetscFunctionReturn(0); 1340 } 1341 1342 /* determine ownership of all rows */ 1343 m = M/size + ((M % size) > rank); 1344 ierr = PetscMalloc((size+2)*sizeof(int),&rowners);CHKERRQ(ierr); 1345 ierr = MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 1346 rowners[0] = 0; 1347 for (i=2; i<=size; i++) { 1348 rowners[i] += rowners[i-1]; 1349 } 1350 rstart = rowners[rank]; 1351 rend = rowners[rank+1]; 1352 1353 /* distribute row lengths to all processors */ 1354 ierr = PetscMalloc(2*(rend-rstart+1)*sizeof(int),&ourlens);CHKERRQ(ierr); 1355 offlens = ourlens + (rend-rstart); 1356 if (!rank) { 1357 ierr = PetscMalloc(M*sizeof(int),&rowlengths);CHKERRQ(ierr); 1358 ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 1359 ierr = PetscMalloc(size*sizeof(int),&sndcounts);CHKERRQ(ierr); 1360 for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i]; 1361 ierr = MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 1362 ierr = PetscFree(sndcounts);CHKERRQ(ierr); 1363 } else { 1364 ierr = MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);CHKERRQ(ierr); 1365 } 1366 1367 if (!rank) { 1368 /* calculate the number of nonzeros on each processor */ 1369 ierr = PetscMalloc(size*sizeof(int),&procsnz);CHKERRQ(ierr); 1370 ierr = PetscMemzero(procsnz,size*sizeof(int));CHKERRQ(ierr); 1371 for (i=0; i<size; i++) { 1372 for (j=rowners[i]; j< rowners[i+1]; j++) { 1373 procsnz[i] += rowlengths[j]; 1374 } 1375 } 1376 ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1377 1378 /* determine max buffer needed and allocate it */ 1379 maxnz = 0; 1380 for (i=0; i<size; i++) { 1381 maxnz = PetscMax(maxnz,procsnz[i]); 1382 } 1383 ierr = PetscMalloc(maxnz*sizeof(int),&cols);CHKERRQ(ierr); 1384 1385 /* read in my part of the matrix column indices */ 1386 nz = procsnz[0]; 1387 ierr = PetscMalloc(nz*sizeof(int),&mycols);CHKERRQ(ierr); 1388 ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr); 1389 1390 /* read in every one elses and ship off */ 1391 for (i=1; i<size; i++) { 1392 nz = procsnz[i]; 1393 ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 1394 ierr = MPI_Send(cols,nz,MPI_INT,i,tag,comm);CHKERRQ(ierr); 1395 } 1396 ierr = PetscFree(cols);CHKERRQ(ierr); 1397 } else { 1398 /* determine buffer space needed for message */ 1399 nz = 0; 1400 for (i=0; i<m; i++) { 1401 nz += ourlens[i]; 1402 } 1403 ierr = PetscMalloc((nz+1)*sizeof(int),&mycols);CHKERRQ(ierr); 1404 1405 /* receive message of column indices*/ 1406 ierr = MPI_Recv(mycols,nz,MPI_INT,0,tag,comm,&status);CHKERRQ(ierr); 1407 ierr = MPI_Get_count(&status,MPI_INT,&maxnz);CHKERRQ(ierr); 1408 if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 1409 } 1410 1411 /* loop over local rows, determining number of off diagonal entries */ 1412 ierr = PetscMemzero(offlens,m*sizeof(int));CHKERRQ(ierr); 1413 jj = 0; 1414 for (i=0; i<m; i++) { 1415 for (j=0; j<ourlens[i]; j++) { 1416 if (mycols[jj] < rstart || mycols[jj] >= rend) offlens[i]++; 1417 jj++; 1418 } 1419 } 1420 1421 /* create our matrix */ 1422 for (i=0; i<m; i++) { 1423 ourlens[i] -= offlens[i]; 1424 } 1425 ierr = MatCreate(comm,m,PETSC_DECIDE,M,N,newmat);CHKERRQ(ierr); 1426 ierr = MatSetType(*newmat,type);CHKERRQ(ierr); 1427 ierr = MatMPIDenseSetPreallocation(*newmat,PETSC_NULL);CHKERRQ(ierr); 1428 A = *newmat; 1429 for (i=0; i<m; i++) { 1430 ourlens[i] += offlens[i]; 1431 } 1432 1433 if (!rank) { 1434 ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 1435 1436 /* read in my part of the matrix numerical values */ 1437 nz = procsnz[0]; 1438 ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1439 1440 /* insert into matrix */ 1441 jj = rstart; 1442 smycols = mycols; 1443 svals = vals; 1444 for (i=0; i<m; i++) { 1445 ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1446 smycols += ourlens[i]; 1447 svals += ourlens[i]; 1448 jj++; 1449 } 1450 1451 /* read in other processors and ship out */ 1452 for (i=1; i<size; i++) { 1453 nz = procsnz[i]; 1454 ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 1455 ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);CHKERRQ(ierr); 1456 } 1457 ierr = PetscFree(procsnz);CHKERRQ(ierr); 1458 } else { 1459 /* receive numeric values */ 1460 ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 1461 1462 /* receive message of values*/ 1463 ierr = MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);CHKERRQ(ierr); 1464 ierr = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr); 1465 if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file"); 1466 1467 /* insert into matrix */ 1468 jj = rstart; 1469 smycols = mycols; 1470 svals = vals; 1471 for (i=0; i<m; i++) { 1472 ierr = MatSetValues(A,1,&jj,ourlens[i],smycols,svals,INSERT_VALUES);CHKERRQ(ierr); 1473 smycols += ourlens[i]; 1474 svals += ourlens[i]; 1475 jj++; 1476 } 1477 } 1478 ierr = PetscFree(ourlens);CHKERRQ(ierr); 1479 ierr = PetscFree(vals);CHKERRQ(ierr); 1480 ierr = PetscFree(mycols);CHKERRQ(ierr); 1481 ierr = PetscFree(rowners);CHKERRQ(ierr); 1482 1483 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1484 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1485 PetscFunctionReturn(0); 1486 } 1487 1488 1489 1490 1491 1492