1 #define PETSCMAT_DLL 2 3 /* 4 This provides a matrix that consists of Mats 5 */ 6 7 #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 8 #include "src/mat/impls/baij/seq/baij.h" /* use the common AIJ data-structure */ 9 10 #define CHUNKSIZE 15 11 12 typedef struct { 13 SEQAIJHEADER(Mat); 14 SEQBAIJHEADER; 15 16 Vec left,right; /* dummy vectors to perform local parts of product */ 17 } Mat_BlockMat; 18 19 #undef __FUNCT__ 20 #define __FUNCT__ "MatSetValues_BlockMat" 21 PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 22 { 23 Mat_BlockMat *a = (Mat_BlockMat*)A->data; 24 PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 25 PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 26 PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol; 27 PetscErrorCode ierr; 28 PetscInt ridx,cidx; 29 PetscTruth roworiented=a->roworiented; 30 MatScalar value; 31 Mat *ap,*aa = a->a; 32 33 PetscFunctionBegin; 34 for (k=0; k<m; k++) { /* loop over added rows */ 35 row = im[k]; 36 brow = row/bs; 37 if (row < 0) continue; 38 #if defined(PETSC_USE_DEBUG) 39 if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 40 #endif 41 rp = aj + ai[brow]; 42 ap = aa + ai[brow]; 43 rmax = imax[brow]; 44 nrow = ailen[brow]; 45 low = 0; 46 high = nrow; 47 for (l=0; l<n; l++) { /* loop over added columns */ 48 if (in[l] < 0) continue; 49 #if defined(PETSC_USE_DEBUG) 50 if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1); 51 #endif 52 col = in[l]; bcol = col/bs; 53 ridx = row % bs; cidx = col % bs; 54 if (roworiented) { 55 value = v[l + k*n]; 56 } else { 57 value = v[k + l*m]; 58 } 59 if (col <= lastcol) low = 0; else high = nrow; 60 lastcol = col; 61 while (high-low > 7) { 62 t = (low+high)/2; 63 if (rp[t] > bcol) high = t; 64 else low = t; 65 } 66 for (i=low; i<high; i++) { 67 if (rp[i] > bcol) break; 68 if (rp[i] == bcol) { 69 goto noinsert1; 70 } 71 } 72 if (nonew == 1) goto noinsert1; 73 if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 74 MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat); 75 N = nrow++ - 1; high++; 76 /* shift up all the later entries in this row */ 77 for (ii=N; ii>=i; ii--) { 78 rp[ii+1] = rp[ii]; 79 ap[ii+1] = ap[ii]; 80 } 81 if (N>=i) ap[i] = 0; 82 rp[i] = bcol; 83 a->nz++; 84 noinsert1:; 85 if (!*(ap+i)) { 86 ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);CHKERRQ(ierr); 87 } 88 ierr = MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);CHKERRQ(ierr); 89 low = i; 90 } 91 ailen[brow] = nrow; 92 } 93 A->same_nonzero = PETSC_FALSE; 94 PetscFunctionReturn(0); 95 } 96 97 #undef __FUNCT__ 98 #define __FUNCT__ "MatLoad_BlockMat" 99 PetscErrorCode MatLoad_BlockMat(PetscViewer viewer, MatType type,Mat *A) 100 { 101 PetscErrorCode ierr; 102 Mat tmpA; 103 PetscInt i,m,n,bs = 1,ncols; 104 const PetscInt *cols; 105 const PetscScalar *values; 106 107 PetscFunctionBegin; 108 ierr = MatLoad_SeqAIJ(viewer,MATSEQAIJ,&tmpA);CHKERRQ(ierr); 109 110 ierr = MatGetLocalSize(tmpA,&m,&n);CHKERRQ(ierr); 111 ierr = PetscOptionsBegin(PETSC_COMM_SELF,PETSC_NULL,"Options for loading BlockMat matrix","Mat");CHKERRQ(ierr); 112 ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);CHKERRQ(ierr); 113 ierr = PetscOptionsEnd();CHKERRQ(ierr); 114 ierr = MatCreateBlockMat(PETSC_COMM_SELF,m,n,bs,A);CHKERRQ(ierr); 115 for (i=0; i<m; i++) { 116 ierr = MatGetRow(tmpA,i,&ncols,&cols,&values);CHKERRQ(ierr); 117 ierr = MatSetValues(*A,1,&i,ncols,cols,values,INSERT_VALUES);CHKERRQ(ierr); 118 } 119 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 120 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 121 PetscFunctionReturn(0); 122 } 123 124 125 #undef __FUNCT__ 126 #define __FUNCT__ "MatDestroy_BlockMat" 127 PetscErrorCode MatDestroy_BlockMat(Mat mat) 128 { 129 PetscErrorCode ierr; 130 Mat_BlockMat *bmat = (Mat_BlockMat*)mat->data; 131 132 PetscFunctionBegin; 133 if (bmat->right) { 134 ierr = VecDestroy(bmat->right);CHKERRQ(ierr); 135 } 136 if (bmat->left) { 137 ierr = VecDestroy(bmat->left);CHKERRQ(ierr); 138 } 139 ierr = PetscFree(bmat);CHKERRQ(ierr); 140 PetscFunctionReturn(0); 141 } 142 143 #undef __FUNCT__ 144 #define __FUNCT__ "MatMult_BlockMat" 145 PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y) 146 { 147 Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 148 PetscErrorCode ierr; 149 PetscScalar *xx,*yy; 150 PetscInt *aj,i,*ii,jrow,m = A->rmap.n/A->rmap.bs,bs = A->rmap.bs,n,j; 151 Mat *aa; 152 153 PetscFunctionBegin; 154 /* 155 Standard CSR multiply except each entry is a Mat 156 */ 157 ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 158 ierr = VecSet(y,0.0);CHKERRQ(ierr); 159 ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 160 aj = bmat->j; 161 aa = bmat->a; 162 ii = bmat->i; 163 for (i=0; i<m; i++) { 164 jrow = ii[i]; 165 ierr = VecPlaceArray(bmat->left,yy + bs*jrow);CHKERRQ(ierr); 166 n = ii[i+1] - jrow; 167 ierr = VecSet(bmat->left,0.0);CHKERRQ(ierr); 168 for (j=0; j<n; j++) { 169 ierr = VecPlaceArray(bmat->right,xx + bs*aj[jrow]);CHKERRQ(ierr); 170 ierr = MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left); 171 ierr = VecResetArray(bmat->right);CHKERRQ(ierr); 172 jrow++; 173 } 174 ierr = VecResetArray(bmat->left);CHKERRQ(ierr); 175 } 176 ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 177 ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 178 PetscFunctionReturn(0); 179 } 180 181 #undef __FUNCT__ 182 #define __FUNCT__ "MatMultAdd_BlockMat" 183 PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 184 { 185 Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 186 PetscErrorCode ierr; 187 188 PetscFunctionBegin; 189 PetscFunctionReturn(0); 190 } 191 192 #undef __FUNCT__ 193 #define __FUNCT__ "MatMultTranspose_BlockMat" 194 PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y) 195 { 196 Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 197 PetscErrorCode ierr; 198 199 PetscFunctionBegin; 200 PetscFunctionReturn(0); 201 } 202 203 #undef __FUNCT__ 204 #define __FUNCT__ "MatMultTransposeAdd_BlockMat" 205 PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z) 206 { 207 Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 208 PetscErrorCode ierr; 209 210 PetscFunctionBegin; 211 PetscFunctionReturn(0); 212 } 213 214 #undef __FUNCT__ 215 #define __FUNCT__ "MatSetBlockSize_BlockMat" 216 PetscErrorCode MatSetBlockSize_BlockMat(Mat A,PetscInt bs) 217 { 218 Mat_BlockMat *bmat = (Mat_BlockMat*)A->data; 219 PetscErrorCode ierr; 220 PetscInt nz = 10,i,m; 221 222 PetscFunctionBegin; 223 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->right);CHKERRQ(ierr); 224 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,PETSC_NULL,&bmat->left);CHKERRQ(ierr); 225 226 A->rmap.bs = A->cmap.bs = bs; 227 228 ierr = PetscMalloc2(A->rmap.n,PetscInt,&bmat->imax,A->rmap.n,PetscInt,&bmat->ilen);CHKERRQ(ierr); 229 for (i=0; i<A->rmap.n; i++) bmat->imax[i] = nz; 230 nz = nz*A->rmap.n; 231 232 233 /* bmat->ilen will count nonzeros in each row so far. */ 234 for (i=0; i<A->rmap.n; i++) { bmat->ilen[i] = 0;} 235 236 /* allocate the matrix space */ 237 ierr = PetscMalloc3(nz,Mat,&bmat->a,nz,PetscInt,&bmat->j,A->rmap.n+1,PetscInt,&bmat->i);CHKERRQ(ierr); 238 bmat->i[0] = 0; 239 for (i=1; i<A->rmap.n+1; i++) { 240 bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1]; 241 } 242 bmat->singlemalloc = PETSC_TRUE; 243 bmat->free_a = PETSC_TRUE; 244 bmat->free_ij = PETSC_TRUE; 245 246 bmat->nz = 0; 247 bmat->maxnz = nz; 248 A->info.nz_unneeded = (double)bmat->maxnz; 249 250 PetscFunctionReturn(0); 251 } 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "MatAssemblyEnd_BlockMat" 255 PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode) 256 { 257 Mat_BlockMat *a = (Mat_BlockMat*)A->data; 258 PetscErrorCode ierr; 259 PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 260 PetscInt m = A->rmap.n/A->rmap.bs,*ip,N,*ailen = a->ilen,rmax = 0; 261 Mat *aa = a->a,*ap; 262 PetscReal ratio=0.6; 263 264 PetscFunctionBegin; 265 if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 266 267 if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 268 for (i=1; i<m; i++) { 269 /* move each row back by the amount of empty slots (fshift) before it*/ 270 fshift += imax[i-1] - ailen[i-1]; 271 rmax = PetscMax(rmax,ailen[i]); 272 if (fshift) { 273 ip = aj + ai[i] ; 274 ap = aa + ai[i] ; 275 N = ailen[i]; 276 for (j=0; j<N; j++) { 277 ip[j-fshift] = ip[j]; 278 ap[j-fshift] = ap[j]; 279 } 280 } 281 ai[i] = ai[i-1] + ailen[i-1]; 282 } 283 if (m) { 284 fshift += imax[m-1] - ailen[m-1]; 285 ai[m] = ai[m-1] + ailen[m-1]; 286 } 287 /* reset ilen and imax for each row */ 288 for (i=0; i<m; i++) { 289 ailen[i] = imax[i] = ai[i+1] - ai[i]; 290 } 291 a->nz = ai[m]; 292 293 ierr = PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap.n,fshift,a->nz);CHKERRQ(ierr); 294 ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr); 295 ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr); 296 a->reallocs = 0; 297 A->info.nz_unneeded = (double)fshift; 298 a->rmax = rmax; 299 300 A->same_nonzero = PETSC_TRUE; 301 PetscFunctionReturn(0); 302 } 303 304 static struct _MatOps MatOps_Values = {MatSetValues_BlockMat, 305 0, 306 0, 307 MatMult_BlockMat, 308 /* 4*/ MatMultAdd_BlockMat, 309 MatMultTranspose_BlockMat, 310 MatMultTransposeAdd_BlockMat, 311 0, 312 0, 313 0, 314 /*10*/ 0, 315 0, 316 0, 317 0, 318 0, 319 /*15*/ 0, 320 0, 321 0, 322 0, 323 0, 324 /*20*/ 0, 325 MatAssemblyEnd_BlockMat, 326 0, 327 0, 328 0, 329 /*25*/ 0, 330 0, 331 0, 332 0, 333 0, 334 /*30*/ 0, 335 0, 336 0, 337 0, 338 0, 339 /*35*/ 0, 340 0, 341 0, 342 0, 343 0, 344 /*40*/ 0, 345 0, 346 0, 347 0, 348 0, 349 /*45*/ 0, 350 0, 351 0, 352 0, 353 0, 354 /*50*/ MatSetBlockSize_BlockMat, 355 0, 356 0, 357 0, 358 0, 359 /*55*/ 0, 360 0, 361 0, 362 0, 363 0, 364 /*60*/ 0, 365 MatDestroy_BlockMat, 366 0, 367 0, 368 0, 369 /*65*/ 0, 370 0, 371 0, 372 0, 373 0, 374 /*70*/ 0, 375 0, 376 0, 377 0, 378 0, 379 /*75*/ 0, 380 0, 381 0, 382 0, 383 0, 384 /*80*/ 0, 385 0, 386 0, 387 0, 388 MatLoad_BlockMat, 389 /*85*/ 0, 390 0, 391 0, 392 0, 393 0, 394 /*90*/ 0, 395 0, 396 0, 397 0, 398 0, 399 /*95*/ 0, 400 0, 401 0, 402 0}; 403 404 /*MC 405 MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix 406 consisting of (usually) sparse blocks. 407 408 Level: advanced 409 410 .seealso: MatCreateBlockMat() 411 412 M*/ 413 414 EXTERN_C_BEGIN 415 #undef __FUNCT__ 416 #define __FUNCT__ "MatCreate_BlockMat" 417 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_BlockMat(Mat A) 418 { 419 Mat_BlockMat *b; 420 PetscErrorCode ierr; 421 422 PetscFunctionBegin; 423 ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 424 ierr = PetscNew(Mat_BlockMat,&b);CHKERRQ(ierr); 425 426 A->data = (void*)b; 427 428 ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr); 429 ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr); 430 431 A->assembled = PETSC_TRUE; 432 A->preallocated = PETSC_FALSE; 433 ierr = PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);CHKERRQ(ierr); 434 PetscFunctionReturn(0); 435 } 436 EXTERN_C_END 437 438 #undef __FUNCT__ 439 #define __FUNCT__ "MatCreateBlockMat" 440 /*@C 441 MatCreateBlockMat - Creates a new matrix based sparse Mat storage 442 443 Collective on MPI_Comm 444 445 Input Parameters: 446 + comm - MPI communicator 447 . m - number of rows 448 . n - number of columns 449 - bs - size of each submatrix 450 451 452 Output Parameter: 453 . A - the matrix 454 455 Level: intermediate 456 457 PETSc requires that matrices and vectors being used for certain 458 operations are partitioned accordingly. For example, when 459 creating a bmat matrix, A, that supports parallel matrix-vector 460 products using MatMult(A,x,y) the user should set the number 461 of local matrix rows to be the number of local elements of the 462 corresponding result vector, y. Note that this is information is 463 required for use of the matrix interface routines, even though 464 the bmat matrix may not actually be physically partitioned. 465 For example, 466 467 .keywords: matrix, bmat, create 468 469 .seealso: MATBLOCKMAT 470 @*/ 471 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,Mat *A) 472 { 473 PetscErrorCode ierr; 474 475 PetscFunctionBegin; 476 ierr = MatCreate(comm,A);CHKERRQ(ierr); 477 ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 478 ierr = MatSetType(*A,MATBLOCKMAT);CHKERRQ(ierr); 479 ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 480 PetscFunctionReturn(0); 481 } 482 483 484 485