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