1d4002b98SHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I "petscmat.h" I*/ 2d4002b98SHong Zhang #include <../src/mat/impls/sell/mpi/mpisell.h> /*I "petscmat.h" I*/ 3d4002b98SHong Zhang #include <petsc/private/vecimpl.h> 4d4002b98SHong Zhang #include <petsc/private/isimpl.h> 5d4002b98SHong Zhang #include <petscblaslapack.h> 6d4002b98SHong Zhang #include <petscsf.h> 7d4002b98SHong Zhang 8d4002b98SHong Zhang /*MC 9d4002b98SHong Zhang MATSELL - MATSELL = "sell" - A matrix type to be used for sparse matrices. 10d4002b98SHong Zhang 11d4002b98SHong Zhang This matrix type is identical to MATSEQSELL when constructed with a single process communicator, 12d4002b98SHong Zhang and MATMPISELL otherwise. As a result, for single process communicators, 13d4002b98SHong Zhang MatSeqSELLSetPreallocation is supported, and similarly MatMPISELLSetPreallocation is supported 14d4002b98SHong Zhang for communicators controlling multiple processes. It is recommended that you call both of 15d4002b98SHong Zhang the above preallocation routines for simplicity. 16d4002b98SHong Zhang 17d4002b98SHong Zhang Options Database Keys: 18d4002b98SHong Zhang . -mat_type sell - sets the matrix type to "sell" during a call to MatSetFromOptions() 19d4002b98SHong Zhang 2095452b02SPatrick Sanan Developer Notes: 2195452b02SPatrick Sanan Subclasses include MATSELLCUSP, MATSELLCUSPARSE, MATSELLPERM, MATSELLCRL, and also automatically switches over to use inodes when 22d4002b98SHong Zhang enough exist. 23d4002b98SHong Zhang 24d4002b98SHong Zhang Level: beginner 25d4002b98SHong Zhang 26d4002b98SHong Zhang .seealso: MatCreateSELL(), MatCreateSeqSELL(), MATSEQSELL, MATMPISELL 27d4002b98SHong Zhang M*/ 28d4002b98SHong Zhang 29d4002b98SHong Zhang PetscErrorCode MatDiagonalSet_MPISELL(Mat Y,Vec D,InsertMode is) 30d4002b98SHong Zhang { 31d4002b98SHong Zhang PetscErrorCode ierr; 32d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)Y->data; 33d4002b98SHong Zhang 34d4002b98SHong Zhang PetscFunctionBegin; 35d4002b98SHong Zhang if (Y->assembled && Y->rmap->rstart == Y->cmap->rstart && Y->rmap->rend == Y->cmap->rend) { 36d4002b98SHong Zhang ierr = MatDiagonalSet(sell->A,D,is);CHKERRQ(ierr); 37d4002b98SHong Zhang } else { 38d4002b98SHong Zhang ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 39d4002b98SHong Zhang } 40d4002b98SHong Zhang PetscFunctionReturn(0); 41d4002b98SHong Zhang } 42d4002b98SHong Zhang 43d4002b98SHong Zhang /* 44d4002b98SHong Zhang Local utility routine that creates a mapping from the global column 45d4002b98SHong Zhang number to the local number in the off-diagonal part of the local 46d4002b98SHong Zhang storage of the matrix. When PETSC_USE_CTABLE is used this is scalable at 47d4002b98SHong Zhang a slightly higher hash table cost; without it it is not scalable (each processor 48d4002b98SHong Zhang has an order N integer array but is fast to acess. 49d4002b98SHong Zhang */ 50d4002b98SHong Zhang PetscErrorCode MatCreateColmap_MPISELL_Private(Mat mat) 51d4002b98SHong Zhang { 52d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 53d4002b98SHong Zhang PetscErrorCode ierr; 54d4002b98SHong Zhang PetscInt n=sell->B->cmap->n,i; 55d4002b98SHong Zhang 56d4002b98SHong Zhang PetscFunctionBegin; 57d4002b98SHong Zhang if (!sell->garray) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPISELL Matrix was assembled but is missing garray"); 58d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE) 59d4002b98SHong Zhang ierr = PetscTableCreate(n,mat->cmap->N+1,&sell->colmap);CHKERRQ(ierr); 60d4002b98SHong Zhang for (i=0; i<n; i++) { 61d4002b98SHong Zhang ierr = PetscTableAdd(sell->colmap,sell->garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr); 62d4002b98SHong Zhang } 63d4002b98SHong Zhang #else 64d4002b98SHong Zhang ierr = PetscCalloc1(mat->cmap->N+1,&sell->colmap);CHKERRQ(ierr); 65d4002b98SHong Zhang ierr = PetscLogObjectMemory((PetscObject)mat,(mat->cmap->N+1)*sizeof(PetscInt));CHKERRQ(ierr); 66d4002b98SHong Zhang for (i=0; i<n; i++) sell->colmap[sell->garray[i]] = i+1; 67d4002b98SHong Zhang #endif 68d4002b98SHong Zhang PetscFunctionReturn(0); 69d4002b98SHong Zhang } 70d4002b98SHong Zhang 71d4002b98SHong Zhang #define MatSetValues_SeqSELL_A_Private(row,col,value,addv,orow,ocol) \ 72d4002b98SHong Zhang { \ 73d4002b98SHong Zhang if (col <= lastcol1) low1 = 0; \ 74d4002b98SHong Zhang else high1 = nrow1; \ 75d4002b98SHong Zhang lastcol1 = col; \ 76d4002b98SHong Zhang while (high1-low1 > 5) { \ 77d4002b98SHong Zhang t = (low1+high1)/2; \ 78d4002b98SHong Zhang if (*(cp1+8*t) > col) high1 = t; \ 79d4002b98SHong Zhang else low1 = t; \ 80d4002b98SHong Zhang } \ 81d4002b98SHong Zhang for (_i=low1; _i<high1; _i++) { \ 82d4002b98SHong Zhang if (*(cp1+8*_i) > col) break; \ 83d4002b98SHong Zhang if (*(cp1+8*_i) == col) { \ 84d4002b98SHong Zhang if (addv == ADD_VALUES) *(vp1+8*_i) += value; \ 85d4002b98SHong Zhang else *(vp1+8*_i) = value; \ 86d4002b98SHong Zhang goto a_noinsert; \ 87d4002b98SHong Zhang } \ 88d4002b98SHong Zhang } \ 89d4002b98SHong Zhang if (value == 0.0 && ignorezeroentries) {low1 = 0; high1 = nrow1;goto a_noinsert;} \ 90d4002b98SHong Zhang if (nonew == 1) {low1 = 0; high1 = nrow1; goto a_noinsert;} \ 91d4002b98SHong Zhang if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at global row/column (%D, %D) into matrix", orow, ocol); \ 92d4002b98SHong Zhang MatSeqXSELLReallocateSELL(A,am,1,nrow1,a->sliidx,row/8,row,col,a->colidx,a->val,cp1,vp1,nonew,MatScalar); \ 93d4002b98SHong Zhang /* shift up all the later entries in this row */ \ 94d4002b98SHong Zhang for (ii=nrow1-1; ii>=_i; ii--) { \ 95d4002b98SHong Zhang *(cp1+8*(ii+1)) = *(cp1+8*ii); \ 96d4002b98SHong Zhang *(vp1+8*(ii+1)) = *(vp1+8*ii); \ 97d4002b98SHong Zhang } \ 98d4002b98SHong Zhang *(cp1+8*_i) = col; \ 99d4002b98SHong Zhang *(vp1+8*_i) = value; \ 100d4002b98SHong Zhang a->nz++; nrow1++; A->nonzerostate++; \ 101d4002b98SHong Zhang a_noinsert: ; \ 102d4002b98SHong Zhang a->rlen[row] = nrow1; \ 103d4002b98SHong Zhang } 104d4002b98SHong Zhang 105d4002b98SHong Zhang #define MatSetValues_SeqSELL_B_Private(row,col,value,addv,orow,ocol) \ 106d4002b98SHong Zhang { \ 107d4002b98SHong Zhang if (col <= lastcol2) low2 = 0; \ 108d4002b98SHong Zhang else high2 = nrow2; \ 109d4002b98SHong Zhang lastcol2 = col; \ 110d4002b98SHong Zhang while (high2-low2 > 5) { \ 111d4002b98SHong Zhang t = (low2+high2)/2; \ 112d4002b98SHong Zhang if (*(cp2+8*t) > col) high2 = t; \ 113d4002b98SHong Zhang else low2 = t; \ 114d4002b98SHong Zhang } \ 115d4002b98SHong Zhang for (_i=low2; _i<high2; _i++) { \ 116d4002b98SHong Zhang if (*(cp2+8*_i) > col) break; \ 117d4002b98SHong Zhang if (*(cp2+8*_i) == col) { \ 118d4002b98SHong Zhang if (addv == ADD_VALUES) *(vp2+8*_i) += value; \ 119d4002b98SHong Zhang else *(vp2+8*_i) = value; \ 120d4002b98SHong Zhang goto b_noinsert; \ 121d4002b98SHong Zhang } \ 122d4002b98SHong Zhang } \ 123d4002b98SHong Zhang if (value == 0.0 && ignorezeroentries) {low2 = 0; high2 = nrow2; goto b_noinsert;} \ 124d4002b98SHong Zhang if (nonew == 1) {low2 = 0; high2 = nrow2; goto b_noinsert;} \ 125d4002b98SHong Zhang if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at global row/column (%D, %D) into matrix", orow, ocol); \ 126d4002b98SHong Zhang MatSeqXSELLReallocateSELL(B,bm,1,nrow2,b->sliidx,row/8,row,col,b->colidx,b->val,cp2,vp2,nonew,MatScalar); \ 127d4002b98SHong Zhang /* shift up all the later entries in this row */ \ 128d4002b98SHong Zhang for (ii=nrow2-1; ii>=_i; ii--) { \ 129d4002b98SHong Zhang *(cp2+8*(ii+1)) = *(cp2+8*ii); \ 130d4002b98SHong Zhang *(vp2+8*(ii+1)) = *(vp2+8*ii); \ 131d4002b98SHong Zhang } \ 132d4002b98SHong Zhang *(cp2+8*_i) = col; \ 133d4002b98SHong Zhang *(vp2+8*_i) = value; \ 134d4002b98SHong Zhang b->nz++; nrow2++; B->nonzerostate++; \ 135d4002b98SHong Zhang b_noinsert: ; \ 136d4002b98SHong Zhang b->rlen[row] = nrow2; \ 137d4002b98SHong Zhang } 138d4002b98SHong Zhang 139d4002b98SHong Zhang PetscErrorCode MatSetValues_MPISELL(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv) 140d4002b98SHong Zhang { 141d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 142d4002b98SHong Zhang PetscScalar value; 143d4002b98SHong Zhang PetscErrorCode ierr; 144d4002b98SHong Zhang PetscInt i,j,rstart=mat->rmap->rstart,rend=mat->rmap->rend,shift1,shift2; 145d4002b98SHong Zhang PetscInt cstart=mat->cmap->rstart,cend=mat->cmap->rend,row,col; 146d4002b98SHong Zhang PetscBool roworiented=sell->roworiented; 147d4002b98SHong Zhang 148d4002b98SHong Zhang /* Some Variables required in the macro */ 149d4002b98SHong Zhang Mat A=sell->A; 150d4002b98SHong Zhang Mat_SeqSELL *a=(Mat_SeqSELL*)A->data; 151d4002b98SHong Zhang PetscBool ignorezeroentries=a->ignorezeroentries,found; 152d4002b98SHong Zhang Mat B=sell->B; 153d4002b98SHong Zhang Mat_SeqSELL *b=(Mat_SeqSELL*)B->data; 154d4002b98SHong Zhang PetscInt *cp1,*cp2,ii,_i,nrow1,nrow2,low1,high1,low2,high2,t,lastcol1,lastcol2; 155d4002b98SHong Zhang MatScalar *vp1,*vp2; 156d4002b98SHong Zhang 157d4002b98SHong Zhang PetscFunctionBegin; 158d4002b98SHong Zhang for (i=0; i<m; i++) { 159d4002b98SHong Zhang if (im[i] < 0) continue; 160d4002b98SHong Zhang #if defined(PETSC_USE_DEBUG) 161d4002b98SHong Zhang if (im[i] >= mat->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",im[i],mat->rmap->N-1); 162d4002b98SHong Zhang #endif 163d4002b98SHong Zhang if (im[i] >= rstart && im[i] < rend) { 164d4002b98SHong Zhang row = im[i] - rstart; 165d4002b98SHong Zhang lastcol1 = -1; 166d4002b98SHong Zhang shift1 = a->sliidx[row>>3]+(row&0x07); /* starting index of the row */ 167d4002b98SHong Zhang cp1 = a->colidx+shift1; 168d4002b98SHong Zhang vp1 = a->val+shift1; 169d4002b98SHong Zhang nrow1 = a->rlen[row]; 170d4002b98SHong Zhang low1 = 0; 171d4002b98SHong Zhang high1 = nrow1; 172d4002b98SHong Zhang lastcol2 = -1; 173d4002b98SHong Zhang shift2 = b->sliidx[row>>3]+(row&0x07); /* starting index of the row */ 174d4002b98SHong Zhang cp2 = b->colidx+shift2; 175d4002b98SHong Zhang vp2 = b->val+shift2; 176d4002b98SHong Zhang nrow2 = b->rlen[row]; 177d4002b98SHong Zhang low2 = 0; 178d4002b98SHong Zhang high2 = nrow2; 179d4002b98SHong Zhang 180d4002b98SHong Zhang for (j=0; j<n; j++) { 181d4002b98SHong Zhang if (roworiented) value = v[i*n+j]; 182d4002b98SHong Zhang else value = v[i+j*m]; 183d4002b98SHong Zhang if (ignorezeroentries && value == 0.0 && (addv == ADD_VALUES)) continue; 184d4002b98SHong Zhang if (in[j] >= cstart && in[j] < cend) { 185d4002b98SHong Zhang col = in[j] - cstart; 186d4002b98SHong Zhang MatSetValue_SeqSELL_Private(A,row,col,value,addv,im[i],in[j],cp1,vp1,lastcol1,low1,high1); /* set one value */ 187d4002b98SHong Zhang } else if (in[j] < 0) continue; 188d4002b98SHong Zhang #if defined(PETSC_USE_DEBUG) 189d4002b98SHong Zhang else if (in[j] >= mat->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[j],mat->cmap->N-1); 190d4002b98SHong Zhang #endif 191d4002b98SHong Zhang else { 192d4002b98SHong Zhang if (mat->was_assembled) { 193d4002b98SHong Zhang if (!sell->colmap) { 194d4002b98SHong Zhang ierr = MatCreateColmap_MPISELL_Private(mat);CHKERRQ(ierr); 195d4002b98SHong Zhang } 196d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE) 197d4002b98SHong Zhang ierr = PetscTableFind(sell->colmap,in[j]+1,&col);CHKERRQ(ierr); 198d4002b98SHong Zhang col--; 199d4002b98SHong Zhang #else 200d4002b98SHong Zhang col = sell->colmap[in[j]] - 1; 201d4002b98SHong Zhang #endif 202d4002b98SHong Zhang if (col < 0 && !((Mat_SeqSELL*)(sell->B->data))->nonew) { 203d4002b98SHong Zhang ierr = MatDisAssemble_MPISELL(mat);CHKERRQ(ierr); 204d4002b98SHong Zhang col = in[j]; 205d4002b98SHong Zhang /* Reinitialize the variables required by MatSetValues_SeqSELL_B_Private() */ 206d4002b98SHong Zhang B = sell->B; 207d4002b98SHong Zhang b = (Mat_SeqSELL*)B->data; 208d4002b98SHong Zhang shift2 = b->sliidx[row>>3]+(row&0x07); /* starting index of the row */ 209d4002b98SHong Zhang cp2 = b->colidx+shift2; 210d4002b98SHong Zhang vp2 = b->val+shift2; 211d4002b98SHong Zhang nrow2 = b->rlen[row]; 212d4002b98SHong Zhang low2 = 0; 213d4002b98SHong Zhang high2 = nrow2; 214d4002b98SHong Zhang } else if (col < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at global row/column (%D, %D) into matrix", im[i], in[j]); 215d4002b98SHong Zhang } else col = in[j]; 216d4002b98SHong Zhang MatSetValue_SeqSELL_Private(B,row,col,value,addv,im[i],in[j],cp2,vp2,lastcol2,low2,high2); /* set one value */ 217d4002b98SHong Zhang } 218d4002b98SHong Zhang } 219d4002b98SHong Zhang } else { 220d4002b98SHong Zhang if (mat->nooffprocentries) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Setting off process row %D even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set",im[i]); 221d4002b98SHong Zhang if (!sell->donotstash) { 222d4002b98SHong Zhang mat->assembled = PETSC_FALSE; 223d4002b98SHong Zhang if (roworiented) { 224d4002b98SHong Zhang ierr = MatStashValuesRow_Private(&mat->stash,im[i],n,in,v+i*n,(PetscBool)(ignorezeroentries && (addv == ADD_VALUES)));CHKERRQ(ierr); 225d4002b98SHong Zhang } else { 226d4002b98SHong Zhang ierr = MatStashValuesCol_Private(&mat->stash,im[i],n,in,v+i,m,(PetscBool)(ignorezeroentries && (addv == ADD_VALUES)));CHKERRQ(ierr); 227d4002b98SHong Zhang } 228d4002b98SHong Zhang } 229d4002b98SHong Zhang } 230d4002b98SHong Zhang } 231d4002b98SHong Zhang PetscFunctionReturn(0); 232d4002b98SHong Zhang } 233d4002b98SHong Zhang 234d4002b98SHong Zhang PetscErrorCode MatGetValues_MPISELL(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[]) 235d4002b98SHong Zhang { 236d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 237d4002b98SHong Zhang PetscErrorCode ierr; 238d4002b98SHong Zhang PetscInt i,j,rstart=mat->rmap->rstart,rend=mat->rmap->rend; 239d4002b98SHong Zhang PetscInt cstart=mat->cmap->rstart,cend=mat->cmap->rend,row,col; 240d4002b98SHong Zhang 241d4002b98SHong Zhang PetscFunctionBegin; 242d4002b98SHong Zhang for (i=0; i<m; i++) { 243d4002b98SHong Zhang if (idxm[i] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",idxm[i]);*/ 244d4002b98SHong Zhang if (idxm[i] >= mat->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm[i],mat->rmap->N-1); 245d4002b98SHong Zhang if (idxm[i] >= rstart && idxm[i] < rend) { 246d4002b98SHong Zhang row = idxm[i] - rstart; 247d4002b98SHong Zhang for (j=0; j<n; j++) { 248d4002b98SHong Zhang if (idxn[j] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",idxn[j]); */ 249d4002b98SHong Zhang if (idxn[j] >= mat->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",idxn[j],mat->cmap->N-1); 250d4002b98SHong Zhang if (idxn[j] >= cstart && idxn[j] < cend) { 251d4002b98SHong Zhang col = idxn[j] - cstart; 252d4002b98SHong Zhang ierr = MatGetValues(sell->A,1,&row,1,&col,v+i*n+j);CHKERRQ(ierr); 253d4002b98SHong Zhang } else { 254d4002b98SHong Zhang if (!sell->colmap) { 255d4002b98SHong Zhang ierr = MatCreateColmap_MPISELL_Private(mat);CHKERRQ(ierr); 256d4002b98SHong Zhang } 257d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE) 258d4002b98SHong Zhang ierr = PetscTableFind(sell->colmap,idxn[j]+1,&col);CHKERRQ(ierr); 259d4002b98SHong Zhang col--; 260d4002b98SHong Zhang #else 261d4002b98SHong Zhang col = sell->colmap[idxn[j]] - 1; 262d4002b98SHong Zhang #endif 263d4002b98SHong Zhang if ((col < 0) || (sell->garray[col] != idxn[j])) *(v+i*n+j) = 0.0; 264d4002b98SHong Zhang else { 265d4002b98SHong Zhang ierr = MatGetValues(sell->B,1,&row,1,&col,v+i*n+j);CHKERRQ(ierr); 266d4002b98SHong Zhang } 267d4002b98SHong Zhang } 268d4002b98SHong Zhang } 269d4002b98SHong Zhang } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported"); 270d4002b98SHong Zhang } 271d4002b98SHong Zhang PetscFunctionReturn(0); 272d4002b98SHong Zhang } 273d4002b98SHong Zhang 274d4002b98SHong Zhang extern PetscErrorCode MatMultDiagonalBlock_MPISELL(Mat,Vec,Vec); 275d4002b98SHong Zhang 276d4002b98SHong Zhang PetscErrorCode MatAssemblyBegin_MPISELL(Mat mat,MatAssemblyType mode) 277d4002b98SHong Zhang { 278d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 279d4002b98SHong Zhang PetscErrorCode ierr; 280d4002b98SHong Zhang PetscInt nstash,reallocs; 281d4002b98SHong Zhang 282d4002b98SHong Zhang PetscFunctionBegin; 283d4002b98SHong Zhang if (sell->donotstash || mat->nooffprocentries) PetscFunctionReturn(0); 284d4002b98SHong Zhang 285d4002b98SHong Zhang ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr); 286d4002b98SHong Zhang ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr); 287d4002b98SHong Zhang ierr = PetscInfo2(sell->A,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr); 288d4002b98SHong Zhang PetscFunctionReturn(0); 289d4002b98SHong Zhang } 290d4002b98SHong Zhang 291d4002b98SHong Zhang PetscErrorCode MatAssemblyEnd_MPISELL(Mat mat,MatAssemblyType mode) 292d4002b98SHong Zhang { 293d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 294d4002b98SHong Zhang PetscErrorCode ierr; 295d4002b98SHong Zhang PetscMPIInt n; 296d4002b98SHong Zhang PetscInt i,flg; 297d4002b98SHong Zhang PetscInt *row,*col; 298d4002b98SHong Zhang PetscScalar *val; 299d4002b98SHong Zhang PetscBool other_disassembled; 300d4002b98SHong Zhang /* do not use 'b = (Mat_SeqSELL*)sell->B->data' as B can be reset in disassembly */ 301d4002b98SHong Zhang PetscFunctionBegin; 302d4002b98SHong Zhang if (!sell->donotstash && !mat->nooffprocentries) { 303d4002b98SHong Zhang while (1) { 304d4002b98SHong Zhang ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 305d4002b98SHong Zhang if (!flg) break; 306d4002b98SHong Zhang 307d4002b98SHong Zhang for (i=0; i<n; i++) { /* assemble one by one */ 308d4002b98SHong Zhang ierr = MatSetValues_MPISELL(mat,1,row+i,1,col+i,val+i,mat->insertmode);CHKERRQ(ierr); 309d4002b98SHong Zhang } 310d4002b98SHong Zhang } 311d4002b98SHong Zhang ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr); 312d4002b98SHong Zhang } 313d4002b98SHong Zhang ierr = MatAssemblyBegin(sell->A,mode);CHKERRQ(ierr); 314d4002b98SHong Zhang ierr = MatAssemblyEnd(sell->A,mode);CHKERRQ(ierr); 315d4002b98SHong Zhang 316d4002b98SHong Zhang /* 317d4002b98SHong Zhang determine if any processor has disassembled, if so we must 318d4002b98SHong Zhang also disassemble ourselfs, in order that we may reassemble. 319d4002b98SHong Zhang */ 320d4002b98SHong Zhang /* 321d4002b98SHong Zhang if nonzero structure of submatrix B cannot change then we know that 322d4002b98SHong Zhang no processor disassembled thus we can skip this stuff 323d4002b98SHong Zhang */ 324d4002b98SHong Zhang if (!((Mat_SeqSELL*)sell->B->data)->nonew) { 325d4002b98SHong Zhang ierr = MPIU_Allreduce(&mat->was_assembled,&other_disassembled,1,MPIU_BOOL,MPI_PROD,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 326d4002b98SHong Zhang if (mat->was_assembled && !other_disassembled) { 327d4002b98SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatDisAssemble not implemented yet\n"); 328d4002b98SHong Zhang ierr = MatDisAssemble_MPISELL(mat);CHKERRQ(ierr); 329d4002b98SHong Zhang } 330d4002b98SHong Zhang } 331d4002b98SHong Zhang if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) { 332d4002b98SHong Zhang ierr = MatSetUpMultiply_MPISELL(mat);CHKERRQ(ierr); 333d4002b98SHong Zhang } 334d4002b98SHong Zhang /* 335d4002b98SHong Zhang ierr = MatSetOption(sell->B,MAT_USE_INODES,PETSC_FALSE);CHKERRQ(ierr); 336d4002b98SHong Zhang */ 337d4002b98SHong Zhang ierr = MatAssemblyBegin(sell->B,mode);CHKERRQ(ierr); 338d4002b98SHong Zhang ierr = MatAssemblyEnd(sell->B,mode);CHKERRQ(ierr); 339d4002b98SHong Zhang ierr = PetscFree2(sell->rowvalues,sell->rowindices);CHKERRQ(ierr); 340d4002b98SHong Zhang sell->rowvalues = 0; 341d4002b98SHong Zhang ierr = VecDestroy(&sell->diag);CHKERRQ(ierr); 342d4002b98SHong Zhang 343d4002b98SHong Zhang /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */ 344d4002b98SHong Zhang if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqSELL*)(sell->A->data))->nonew) { 345d4002b98SHong Zhang PetscObjectState state = sell->A->nonzerostate + sell->B->nonzerostate; 346d4002b98SHong Zhang ierr = MPIU_Allreduce(&state,&mat->nonzerostate,1,MPIU_INT64,MPI_SUM,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 347d4002b98SHong Zhang } 348d4002b98SHong Zhang PetscFunctionReturn(0); 349d4002b98SHong Zhang } 350d4002b98SHong Zhang 351d4002b98SHong Zhang PetscErrorCode MatZeroEntries_MPISELL(Mat A) 352d4002b98SHong Zhang { 353d4002b98SHong Zhang Mat_MPISELL *l=(Mat_MPISELL*)A->data; 354d4002b98SHong Zhang PetscErrorCode ierr; 355d4002b98SHong Zhang 356d4002b98SHong Zhang PetscFunctionBegin; 357d4002b98SHong Zhang ierr = MatZeroEntries(l->A);CHKERRQ(ierr); 358d4002b98SHong Zhang ierr = MatZeroEntries(l->B);CHKERRQ(ierr); 359d4002b98SHong Zhang PetscFunctionReturn(0); 360d4002b98SHong Zhang } 361d4002b98SHong Zhang 362d4002b98SHong Zhang PetscErrorCode MatMult_MPISELL(Mat A,Vec xx,Vec yy) 363d4002b98SHong Zhang { 364d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 365d4002b98SHong Zhang PetscErrorCode ierr; 366d4002b98SHong Zhang PetscInt nt; 367d4002b98SHong Zhang 368d4002b98SHong Zhang PetscFunctionBegin; 369d4002b98SHong Zhang ierr = VecGetLocalSize(xx,&nt);CHKERRQ(ierr); 370d4002b98SHong Zhang if (nt != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A (%D) and xx (%D)",A->cmap->n,nt); 371d4002b98SHong Zhang ierr = VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 372d4002b98SHong Zhang ierr = (*a->A->ops->mult)(a->A,xx,yy);CHKERRQ(ierr); 373d4002b98SHong Zhang ierr = VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 374d4002b98SHong Zhang ierr = (*a->B->ops->multadd)(a->B,a->lvec,yy,yy);CHKERRQ(ierr); 375d4002b98SHong Zhang PetscFunctionReturn(0); 376d4002b98SHong Zhang } 377d4002b98SHong Zhang 378d4002b98SHong Zhang PetscErrorCode MatMultDiagonalBlock_MPISELL(Mat A,Vec bb,Vec xx) 379d4002b98SHong Zhang { 380d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 381d4002b98SHong Zhang PetscErrorCode ierr; 382d4002b98SHong Zhang 383d4002b98SHong Zhang PetscFunctionBegin; 384d4002b98SHong Zhang ierr = MatMultDiagonalBlock(a->A,bb,xx);CHKERRQ(ierr); 385d4002b98SHong Zhang PetscFunctionReturn(0); 386d4002b98SHong Zhang } 387d4002b98SHong Zhang 388d4002b98SHong Zhang PetscErrorCode MatMultAdd_MPISELL(Mat A,Vec xx,Vec yy,Vec zz) 389d4002b98SHong Zhang { 390d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 391d4002b98SHong Zhang PetscErrorCode ierr; 392d4002b98SHong Zhang 393d4002b98SHong Zhang PetscFunctionBegin; 394d4002b98SHong Zhang ierr = VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 395d4002b98SHong Zhang ierr = (*a->A->ops->multadd)(a->A,xx,yy,zz);CHKERRQ(ierr); 396d4002b98SHong Zhang ierr = VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 397d4002b98SHong Zhang ierr = (*a->B->ops->multadd)(a->B,a->lvec,zz,zz);CHKERRQ(ierr); 398d4002b98SHong Zhang PetscFunctionReturn(0); 399d4002b98SHong Zhang } 400d4002b98SHong Zhang 401d4002b98SHong Zhang PetscErrorCode MatMultTranspose_MPISELL(Mat A,Vec xx,Vec yy) 402d4002b98SHong Zhang { 403d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 404d4002b98SHong Zhang PetscErrorCode ierr; 405d4002b98SHong Zhang 406d4002b98SHong Zhang PetscFunctionBegin; 407d4002b98SHong Zhang /* do nondiagonal part */ 408d4002b98SHong Zhang ierr = (*a->B->ops->multtranspose)(a->B,xx,a->lvec);CHKERRQ(ierr); 409d4002b98SHong Zhang /* do local part */ 410d4002b98SHong Zhang ierr = (*a->A->ops->multtranspose)(a->A,xx,yy);CHKERRQ(ierr); 411*a29b4eb7SJunchao Zhang /* add partial results together */ 412d4002b98SHong Zhang ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 413d4002b98SHong Zhang ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 414d4002b98SHong Zhang PetscFunctionReturn(0); 415d4002b98SHong Zhang } 416d4002b98SHong Zhang 417d4002b98SHong Zhang PetscErrorCode MatIsTranspose_MPISELL(Mat Amat,Mat Bmat,PetscReal tol,PetscBool *f) 418d4002b98SHong Zhang { 419d4002b98SHong Zhang MPI_Comm comm; 420d4002b98SHong Zhang Mat_MPISELL *Asell=(Mat_MPISELL*)Amat->data,*Bsell; 421d4002b98SHong Zhang Mat Adia=Asell->A,Bdia,Aoff,Boff,*Aoffs,*Boffs; 422d4002b98SHong Zhang IS Me,Notme; 423d4002b98SHong Zhang PetscErrorCode ierr; 424d4002b98SHong Zhang PetscInt M,N,first,last,*notme,i; 425d4002b98SHong Zhang PetscMPIInt size; 426d4002b98SHong Zhang 427d4002b98SHong Zhang PetscFunctionBegin; 428d4002b98SHong Zhang /* Easy test: symmetric diagonal block */ 429d4002b98SHong Zhang Bsell = (Mat_MPISELL*)Bmat->data; Bdia = Bsell->A; 430d4002b98SHong Zhang ierr = MatIsTranspose(Adia,Bdia,tol,f);CHKERRQ(ierr); 431d4002b98SHong Zhang if (!*f) PetscFunctionReturn(0); 432d4002b98SHong Zhang ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 433d4002b98SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 434d4002b98SHong Zhang if (size == 1) PetscFunctionReturn(0); 435d4002b98SHong Zhang 436d4002b98SHong Zhang /* Hard test: off-diagonal block. This takes a MatCreateSubMatrix. */ 437d4002b98SHong Zhang ierr = MatGetSize(Amat,&M,&N);CHKERRQ(ierr); 438d4002b98SHong Zhang ierr = MatGetOwnershipRange(Amat,&first,&last);CHKERRQ(ierr); 439d4002b98SHong Zhang ierr = PetscMalloc1(N-last+first,¬me);CHKERRQ(ierr); 440d4002b98SHong Zhang for (i=0; i<first; i++) notme[i] = i; 441d4002b98SHong Zhang for (i=last; i<M; i++) notme[i-last+first] = i; 442d4002b98SHong Zhang ierr = ISCreateGeneral(MPI_COMM_SELF,N-last+first,notme,PETSC_COPY_VALUES,&Notme);CHKERRQ(ierr); 443d4002b98SHong Zhang ierr = ISCreateStride(MPI_COMM_SELF,last-first,first,1,&Me);CHKERRQ(ierr); 444d4002b98SHong Zhang ierr = MatCreateSubMatrices(Amat,1,&Me,&Notme,MAT_INITIAL_MATRIX,&Aoffs);CHKERRQ(ierr); 445d4002b98SHong Zhang Aoff = Aoffs[0]; 446d4002b98SHong Zhang ierr = MatCreateSubMatrices(Bmat,1,&Notme,&Me,MAT_INITIAL_MATRIX,&Boffs);CHKERRQ(ierr); 447d4002b98SHong Zhang Boff = Boffs[0]; 448d4002b98SHong Zhang ierr = MatIsTranspose(Aoff,Boff,tol,f);CHKERRQ(ierr); 449d4002b98SHong Zhang ierr = MatDestroyMatrices(1,&Aoffs);CHKERRQ(ierr); 450d4002b98SHong Zhang ierr = MatDestroyMatrices(1,&Boffs);CHKERRQ(ierr); 451d4002b98SHong Zhang ierr = ISDestroy(&Me);CHKERRQ(ierr); 452d4002b98SHong Zhang ierr = ISDestroy(&Notme);CHKERRQ(ierr); 453d4002b98SHong Zhang ierr = PetscFree(notme);CHKERRQ(ierr); 454d4002b98SHong Zhang PetscFunctionReturn(0); 455d4002b98SHong Zhang } 456d4002b98SHong Zhang 457d4002b98SHong Zhang PetscErrorCode MatMultTransposeAdd_MPISELL(Mat A,Vec xx,Vec yy,Vec zz) 458d4002b98SHong Zhang { 459d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 460d4002b98SHong Zhang PetscErrorCode ierr; 461d4002b98SHong Zhang 462d4002b98SHong Zhang PetscFunctionBegin; 463d4002b98SHong Zhang /* do nondiagonal part */ 464d4002b98SHong Zhang ierr = (*a->B->ops->multtranspose)(a->B,xx,a->lvec);CHKERRQ(ierr); 465d4002b98SHong Zhang /* do local part */ 466d4002b98SHong Zhang ierr = (*a->A->ops->multtransposeadd)(a->A,xx,yy,zz);CHKERRQ(ierr); 467e4a140f6SJunchao Zhang /* add partial results together */ 468e4a140f6SJunchao Zhang ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 469d4002b98SHong Zhang ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 470d4002b98SHong Zhang PetscFunctionReturn(0); 471d4002b98SHong Zhang } 472d4002b98SHong Zhang 473d4002b98SHong Zhang /* 474d4002b98SHong Zhang This only works correctly for square matrices where the subblock A->A is the 475d4002b98SHong Zhang diagonal block 476d4002b98SHong Zhang */ 477d4002b98SHong Zhang PetscErrorCode MatGetDiagonal_MPISELL(Mat A,Vec v) 478d4002b98SHong Zhang { 479d4002b98SHong Zhang PetscErrorCode ierr; 480d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 481d4002b98SHong Zhang 482d4002b98SHong Zhang PetscFunctionBegin; 483d4002b98SHong Zhang if (A->rmap->N != A->cmap->N) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block"); 484d4002b98SHong Zhang if (A->rmap->rstart != A->cmap->rstart || A->rmap->rend != A->cmap->rend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"row partition must equal col partition"); 485d4002b98SHong Zhang ierr = MatGetDiagonal(a->A,v);CHKERRQ(ierr); 486d4002b98SHong Zhang PetscFunctionReturn(0); 487d4002b98SHong Zhang } 488d4002b98SHong Zhang 489d4002b98SHong Zhang PetscErrorCode MatScale_MPISELL(Mat A,PetscScalar aa) 490d4002b98SHong Zhang { 491d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 492d4002b98SHong Zhang PetscErrorCode ierr; 493d4002b98SHong Zhang 494d4002b98SHong Zhang PetscFunctionBegin; 495d4002b98SHong Zhang ierr = MatScale(a->A,aa);CHKERRQ(ierr); 496d4002b98SHong Zhang ierr = MatScale(a->B,aa);CHKERRQ(ierr); 497d4002b98SHong Zhang PetscFunctionReturn(0); 498d4002b98SHong Zhang } 499d4002b98SHong Zhang 500d4002b98SHong Zhang PetscErrorCode MatDestroy_MPISELL(Mat mat) 501d4002b98SHong Zhang { 502d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 503d4002b98SHong Zhang PetscErrorCode ierr; 504d4002b98SHong Zhang 505d4002b98SHong Zhang PetscFunctionBegin; 506d4002b98SHong Zhang #if defined(PETSC_USE_LOG) 507d4002b98SHong Zhang PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D",mat->rmap->N,mat->cmap->N); 508d4002b98SHong Zhang #endif 509d4002b98SHong Zhang ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr); 510d4002b98SHong Zhang ierr = VecDestroy(&sell->diag);CHKERRQ(ierr); 511d4002b98SHong Zhang ierr = MatDestroy(&sell->A);CHKERRQ(ierr); 512d4002b98SHong Zhang ierr = MatDestroy(&sell->B);CHKERRQ(ierr); 513d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE) 514d4002b98SHong Zhang ierr = PetscTableDestroy(&sell->colmap);CHKERRQ(ierr); 515d4002b98SHong Zhang #else 516d4002b98SHong Zhang ierr = PetscFree(sell->colmap);CHKERRQ(ierr); 517d4002b98SHong Zhang #endif 518d4002b98SHong Zhang ierr = PetscFree(sell->garray);CHKERRQ(ierr); 519d4002b98SHong Zhang ierr = VecDestroy(&sell->lvec);CHKERRQ(ierr); 520d4002b98SHong Zhang ierr = VecScatterDestroy(&sell->Mvctx);CHKERRQ(ierr); 521d4002b98SHong Zhang ierr = PetscFree2(sell->rowvalues,sell->rowindices);CHKERRQ(ierr); 522d4002b98SHong Zhang ierr = PetscFree(sell->ld);CHKERRQ(ierr); 523d4002b98SHong Zhang ierr = PetscFree(mat->data);CHKERRQ(ierr); 524d4002b98SHong Zhang 525d4002b98SHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 526d4002b98SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C",NULL);CHKERRQ(ierr); 527d4002b98SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C",NULL);CHKERRQ(ierr); 528d4002b98SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatIsTranspose_C",NULL);CHKERRQ(ierr); 529d4002b98SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPISELLSetPreallocation_C",NULL);CHKERRQ(ierr); 530d4002b98SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisell_mpiaij_C",NULL);CHKERRQ(ierr); 531d4002b98SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDiagonalScaleLocal_C",NULL);CHKERRQ(ierr); 532d4002b98SHong Zhang PetscFunctionReturn(0); 533d4002b98SHong Zhang } 534d4002b98SHong Zhang 535d4002b98SHong Zhang #include <petscdraw.h> 536d4002b98SHong Zhang PetscErrorCode MatView_MPISELL_ASCIIorDraworSocket(Mat mat,PetscViewer viewer) 537d4002b98SHong Zhang { 538d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 539d4002b98SHong Zhang PetscErrorCode ierr; 540d4002b98SHong Zhang PetscMPIInt rank=sell->rank,size=sell->size; 541d4002b98SHong Zhang PetscBool isdraw,iascii,isbinary; 542d4002b98SHong Zhang PetscViewer sviewer; 543d4002b98SHong Zhang PetscViewerFormat format; 544d4002b98SHong Zhang 545d4002b98SHong Zhang PetscFunctionBegin; 546d4002b98SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 547d4002b98SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 548d4002b98SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 549d4002b98SHong Zhang if (iascii) { 550d4002b98SHong Zhang ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 551d4002b98SHong Zhang if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 552d4002b98SHong Zhang MatInfo info; 553d4002b98SHong Zhang PetscBool inodes; 554d4002b98SHong Zhang 555d4002b98SHong Zhang ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);CHKERRQ(ierr); 556d4002b98SHong Zhang ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 557d4002b98SHong Zhang ierr = MatInodeGetInodeSizes(sell->A,NULL,(PetscInt**)&inodes,NULL);CHKERRQ(ierr); 558d4002b98SHong Zhang ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 559d4002b98SHong Zhang if (!inodes) { 560d4002b98SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %D nz %D nz alloced %D mem %D, not using I-node routines\n", 561d4002b98SHong Zhang rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr); 562d4002b98SHong Zhang } else { 563d4002b98SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %D nz %D nz alloced %D mem %D, using I-node routines\n", 564d4002b98SHong Zhang rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,(PetscInt)info.memory);CHKERRQ(ierr); 565d4002b98SHong Zhang } 566d4002b98SHong Zhang ierr = MatGetInfo(sell->A,MAT_LOCAL,&info);CHKERRQ(ierr); 567d4002b98SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);CHKERRQ(ierr); 568d4002b98SHong Zhang ierr = MatGetInfo(sell->B,MAT_LOCAL,&info);CHKERRQ(ierr); 569d4002b98SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);CHKERRQ(ierr); 570d4002b98SHong Zhang ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 571d4002b98SHong Zhang ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 572d4002b98SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n");CHKERRQ(ierr); 573d4002b98SHong Zhang ierr = VecScatterView(sell->Mvctx,viewer);CHKERRQ(ierr); 574d4002b98SHong Zhang PetscFunctionReturn(0); 575d4002b98SHong Zhang } else if (format == PETSC_VIEWER_ASCII_INFO) { 576d4002b98SHong Zhang PetscInt inodecount,inodelimit,*inodes; 577d4002b98SHong Zhang ierr = MatInodeGetInodeSizes(sell->A,&inodecount,&inodes,&inodelimit);CHKERRQ(ierr); 578d4002b98SHong Zhang if (inodes) { 579d4002b98SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"using I-node (on process 0) routines: found %D nodes, limit used is %D\n",inodecount,inodelimit);CHKERRQ(ierr); 580d4002b98SHong Zhang } else { 581d4002b98SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"not using I-node (on process 0) routines\n");CHKERRQ(ierr); 582d4002b98SHong Zhang } 583d4002b98SHong Zhang PetscFunctionReturn(0); 584d4002b98SHong Zhang } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 585d4002b98SHong Zhang PetscFunctionReturn(0); 586d4002b98SHong Zhang } 587d4002b98SHong Zhang } else if (isbinary) { 588d4002b98SHong Zhang if (size == 1) { 589d4002b98SHong Zhang ierr = PetscObjectSetName((PetscObject)sell->A,((PetscObject)mat)->name);CHKERRQ(ierr); 590d4002b98SHong Zhang ierr = MatView(sell->A,viewer);CHKERRQ(ierr); 591d4002b98SHong Zhang } else { 592d4002b98SHong Zhang /* ierr = MatView_MPISELL_Binary(mat,viewer);CHKERRQ(ierr); */ 593d4002b98SHong Zhang } 594d4002b98SHong Zhang PetscFunctionReturn(0); 595d4002b98SHong Zhang } else if (isdraw) { 596d4002b98SHong Zhang PetscDraw draw; 597d4002b98SHong Zhang PetscBool isnull; 598d4002b98SHong Zhang ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 599d4002b98SHong Zhang ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 600d4002b98SHong Zhang if (isnull) PetscFunctionReturn(0); 601d4002b98SHong Zhang } 602d4002b98SHong Zhang 603d4002b98SHong Zhang { 604d4002b98SHong Zhang /* assemble the entire matrix onto first processor. */ 605d4002b98SHong Zhang Mat A; 606d4002b98SHong Zhang Mat_SeqSELL *Aloc; 607d4002b98SHong Zhang PetscInt M=mat->rmap->N,N=mat->cmap->N,*acolidx,row,col,i,j; 608d4002b98SHong Zhang MatScalar *aval; 609d4002b98SHong Zhang PetscBool isnonzero; 610d4002b98SHong Zhang 611d4002b98SHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)mat),&A);CHKERRQ(ierr); 612d4002b98SHong Zhang if (!rank) { 613d4002b98SHong Zhang ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr); 614d4002b98SHong Zhang } else { 615d4002b98SHong Zhang ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr); 616d4002b98SHong Zhang } 617d4002b98SHong Zhang /* This is just a temporary matrix, so explicitly using MATMPISELL is probably best */ 618d4002b98SHong Zhang ierr = MatSetType(A,MATMPISELL);CHKERRQ(ierr); 619d4002b98SHong Zhang ierr = MatMPISELLSetPreallocation(A,0,NULL,0,NULL);CHKERRQ(ierr); 620d4002b98SHong Zhang ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 621d4002b98SHong Zhang ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)A);CHKERRQ(ierr); 622d4002b98SHong Zhang 623d4002b98SHong Zhang /* copy over the A part */ 624d4002b98SHong Zhang Aloc = (Mat_SeqSELL*)sell->A->data; 625d4002b98SHong Zhang acolidx = Aloc->colidx; aval = Aloc->val; 626d4002b98SHong Zhang for (i=0; i<Aloc->totalslices; i++) { /* loop over slices */ 627d4002b98SHong Zhang for (j=Aloc->sliidx[i]; j<Aloc->sliidx[i+1]; j++) { 628d4002b98SHong Zhang isnonzero = (PetscBool)((j-Aloc->sliidx[i])/8 < Aloc->rlen[(i<<3)+(j&0x07)]); 629d4002b98SHong Zhang if (isnonzero) { /* check the mask bit */ 630d4002b98SHong Zhang row = (i<<3)+(j&0x07) + mat->rmap->rstart; /* i<<3 is the starting row of this slice */ 631d4002b98SHong Zhang col = *acolidx + mat->rmap->rstart; 632d4002b98SHong Zhang ierr = MatSetValues(A,1,&row,1,&col,aval,INSERT_VALUES);CHKERRQ(ierr); 633d4002b98SHong Zhang } 634d4002b98SHong Zhang aval++; acolidx++; 635d4002b98SHong Zhang } 636d4002b98SHong Zhang } 637d4002b98SHong Zhang 638d4002b98SHong Zhang /* copy over the B part */ 639d4002b98SHong Zhang Aloc = (Mat_SeqSELL*)sell->B->data; 640d4002b98SHong Zhang acolidx = Aloc->colidx; aval = Aloc->val; 641d4002b98SHong Zhang for (i=0; i<Aloc->totalslices; i++) { 642d4002b98SHong Zhang for (j=Aloc->sliidx[i]; j<Aloc->sliidx[i+1]; j++) { 643d4002b98SHong Zhang isnonzero = (PetscBool)((j-Aloc->sliidx[i])/8 < Aloc->rlen[(i<<3)+(j&0x07)]); 644d4002b98SHong Zhang if (isnonzero) { 645d4002b98SHong Zhang row = (i<<3)+(j&0x07) + mat->rmap->rstart; 646d4002b98SHong Zhang col = sell->garray[*acolidx]; 647d4002b98SHong Zhang ierr = MatSetValues(A,1,&row,1,&col,aval,INSERT_VALUES);CHKERRQ(ierr); 648d4002b98SHong Zhang } 649d4002b98SHong Zhang aval++; acolidx++; 650d4002b98SHong Zhang } 651d4002b98SHong Zhang } 652d4002b98SHong Zhang 653d4002b98SHong Zhang ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 654d4002b98SHong Zhang ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 655d4002b98SHong Zhang /* 656d4002b98SHong Zhang Everyone has to call to draw the matrix since the graphics waits are 657d4002b98SHong Zhang synchronized across all processors that share the PetscDraw object 658d4002b98SHong Zhang */ 659d4002b98SHong Zhang ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 660d4002b98SHong Zhang if (!rank) { 661d4002b98SHong Zhang ierr = PetscObjectSetName((PetscObject)((Mat_MPISELL*)(A->data))->A,((PetscObject)mat)->name);CHKERRQ(ierr); 662d4002b98SHong Zhang ierr = MatView_SeqSELL(((Mat_MPISELL*)(A->data))->A,sviewer);CHKERRQ(ierr); 663d4002b98SHong Zhang } 664d4002b98SHong Zhang ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 665d4002b98SHong Zhang ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 666d4002b98SHong Zhang ierr = MatDestroy(&A);CHKERRQ(ierr); 667d4002b98SHong Zhang } 668d4002b98SHong Zhang PetscFunctionReturn(0); 669d4002b98SHong Zhang } 670d4002b98SHong Zhang 671d4002b98SHong Zhang PetscErrorCode MatView_MPISELL(Mat mat,PetscViewer viewer) 672d4002b98SHong Zhang { 673d4002b98SHong Zhang PetscErrorCode ierr; 674d4002b98SHong Zhang PetscBool iascii,isdraw,issocket,isbinary; 675d4002b98SHong Zhang 676d4002b98SHong Zhang PetscFunctionBegin; 677d4002b98SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 678d4002b98SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 679d4002b98SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 680d4002b98SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);CHKERRQ(ierr); 681d4002b98SHong Zhang if (iascii || isdraw || isbinary || issocket) { 682d4002b98SHong Zhang ierr = MatView_MPISELL_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr); 683d4002b98SHong Zhang } 684d4002b98SHong Zhang PetscFunctionReturn(0); 685d4002b98SHong Zhang } 686d4002b98SHong Zhang 687d4002b98SHong Zhang PetscErrorCode MatGetGhosts_MPISELL(Mat mat,PetscInt *nghosts,const PetscInt *ghosts[]) 688d4002b98SHong Zhang { 689d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 690d4002b98SHong Zhang PetscErrorCode ierr; 691d4002b98SHong Zhang 692d4002b98SHong Zhang PetscFunctionBegin; 693d4002b98SHong Zhang ierr = MatGetSize(sell->B,NULL,nghosts);CHKERRQ(ierr); 694d4002b98SHong Zhang if (ghosts) *ghosts = sell->garray; 695d4002b98SHong Zhang PetscFunctionReturn(0); 696d4002b98SHong Zhang } 697d4002b98SHong Zhang 698d4002b98SHong Zhang PetscErrorCode MatGetInfo_MPISELL(Mat matin,MatInfoType flag,MatInfo *info) 699d4002b98SHong Zhang { 700d4002b98SHong Zhang Mat_MPISELL *mat=(Mat_MPISELL*)matin->data; 701d4002b98SHong Zhang Mat A=mat->A,B=mat->B; 702d4002b98SHong Zhang PetscErrorCode ierr; 703d4002b98SHong Zhang PetscReal isend[5],irecv[5]; 704d4002b98SHong Zhang 705d4002b98SHong Zhang PetscFunctionBegin; 706d4002b98SHong Zhang info->block_size = 1.0; 707d4002b98SHong Zhang ierr = MatGetInfo(A,MAT_LOCAL,info);CHKERRQ(ierr); 708d4002b98SHong Zhang 709d4002b98SHong Zhang isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded; 710d4002b98SHong Zhang isend[3] = info->memory; isend[4] = info->mallocs; 711d4002b98SHong Zhang 712d4002b98SHong Zhang ierr = MatGetInfo(B,MAT_LOCAL,info);CHKERRQ(ierr); 713d4002b98SHong Zhang 714d4002b98SHong Zhang isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded; 715d4002b98SHong Zhang isend[3] += info->memory; isend[4] += info->mallocs; 716d4002b98SHong Zhang if (flag == MAT_LOCAL) { 717d4002b98SHong Zhang info->nz_used = isend[0]; 718d4002b98SHong Zhang info->nz_allocated = isend[1]; 719d4002b98SHong Zhang info->nz_unneeded = isend[2]; 720d4002b98SHong Zhang info->memory = isend[3]; 721d4002b98SHong Zhang info->mallocs = isend[4]; 722d4002b98SHong Zhang } else if (flag == MAT_GLOBAL_MAX) { 723d4002b98SHong Zhang ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)matin));CHKERRQ(ierr); 724d4002b98SHong Zhang 725d4002b98SHong Zhang info->nz_used = irecv[0]; 726d4002b98SHong Zhang info->nz_allocated = irecv[1]; 727d4002b98SHong Zhang info->nz_unneeded = irecv[2]; 728d4002b98SHong Zhang info->memory = irecv[3]; 729d4002b98SHong Zhang info->mallocs = irecv[4]; 730d4002b98SHong Zhang } else if (flag == MAT_GLOBAL_SUM) { 731d4002b98SHong Zhang ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)matin));CHKERRQ(ierr); 732d4002b98SHong Zhang 733d4002b98SHong Zhang info->nz_used = irecv[0]; 734d4002b98SHong Zhang info->nz_allocated = irecv[1]; 735d4002b98SHong Zhang info->nz_unneeded = irecv[2]; 736d4002b98SHong Zhang info->memory = irecv[3]; 737d4002b98SHong Zhang info->mallocs = irecv[4]; 738d4002b98SHong Zhang } 739d4002b98SHong Zhang info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */ 740d4002b98SHong Zhang info->fill_ratio_needed = 0; 741d4002b98SHong Zhang info->factor_mallocs = 0; 742d4002b98SHong Zhang PetscFunctionReturn(0); 743d4002b98SHong Zhang } 744d4002b98SHong Zhang 745d4002b98SHong Zhang PetscErrorCode MatSetOption_MPISELL(Mat A,MatOption op,PetscBool flg) 746d4002b98SHong Zhang { 747d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 748d4002b98SHong Zhang PetscErrorCode ierr; 749d4002b98SHong Zhang 750d4002b98SHong Zhang PetscFunctionBegin; 751d4002b98SHong Zhang switch (op) { 752d4002b98SHong Zhang case MAT_NEW_NONZERO_LOCATIONS: 753d4002b98SHong Zhang case MAT_NEW_NONZERO_ALLOCATION_ERR: 754d4002b98SHong Zhang case MAT_UNUSED_NONZERO_LOCATION_ERR: 755d4002b98SHong Zhang case MAT_KEEP_NONZERO_PATTERN: 756d4002b98SHong Zhang case MAT_NEW_NONZERO_LOCATION_ERR: 757d4002b98SHong Zhang case MAT_USE_INODES: 758d4002b98SHong Zhang case MAT_IGNORE_ZERO_ENTRIES: 759d4002b98SHong Zhang MatCheckPreallocated(A,1); 760d4002b98SHong Zhang ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 761d4002b98SHong Zhang ierr = MatSetOption(a->B,op,flg);CHKERRQ(ierr); 762d4002b98SHong Zhang break; 763d4002b98SHong Zhang case MAT_ROW_ORIENTED: 764d4002b98SHong Zhang MatCheckPreallocated(A,1); 765d4002b98SHong Zhang a->roworiented = flg; 766d4002b98SHong Zhang 767d4002b98SHong Zhang ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 768d4002b98SHong Zhang ierr = MatSetOption(a->B,op,flg);CHKERRQ(ierr); 769d4002b98SHong Zhang break; 770d4002b98SHong Zhang case MAT_NEW_DIAGONALS: 771d4002b98SHong Zhang ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 772d4002b98SHong Zhang break; 773d4002b98SHong Zhang case MAT_IGNORE_OFF_PROC_ENTRIES: 774d4002b98SHong Zhang a->donotstash = flg; 775d4002b98SHong Zhang break; 776d4002b98SHong Zhang case MAT_SPD: 777d4002b98SHong Zhang A->spd_set = PETSC_TRUE; 778d4002b98SHong Zhang A->spd = flg; 779d4002b98SHong Zhang if (flg) { 780d4002b98SHong Zhang A->symmetric = PETSC_TRUE; 781d4002b98SHong Zhang A->structurally_symmetric = PETSC_TRUE; 782d4002b98SHong Zhang A->symmetric_set = PETSC_TRUE; 783d4002b98SHong Zhang A->structurally_symmetric_set = PETSC_TRUE; 784d4002b98SHong Zhang } 785d4002b98SHong Zhang break; 786d4002b98SHong Zhang case MAT_SYMMETRIC: 787d4002b98SHong Zhang MatCheckPreallocated(A,1); 788d4002b98SHong Zhang ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 789d4002b98SHong Zhang break; 790d4002b98SHong Zhang case MAT_STRUCTURALLY_SYMMETRIC: 791d4002b98SHong Zhang MatCheckPreallocated(A,1); 792d4002b98SHong Zhang ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 793d4002b98SHong Zhang break; 794d4002b98SHong Zhang case MAT_HERMITIAN: 795d4002b98SHong Zhang MatCheckPreallocated(A,1); 796d4002b98SHong Zhang ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 797d4002b98SHong Zhang break; 798d4002b98SHong Zhang case MAT_SYMMETRY_ETERNAL: 799d4002b98SHong Zhang MatCheckPreallocated(A,1); 800d4002b98SHong Zhang ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 801d4002b98SHong Zhang break; 802d4002b98SHong Zhang default: 803d4002b98SHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op); 804d4002b98SHong Zhang } 805d4002b98SHong Zhang PetscFunctionReturn(0); 806d4002b98SHong Zhang } 807d4002b98SHong Zhang 808d4002b98SHong Zhang 809d4002b98SHong Zhang PetscErrorCode MatDiagonalScale_MPISELL(Mat mat,Vec ll,Vec rr) 810d4002b98SHong Zhang { 811d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 812d4002b98SHong Zhang Mat a=sell->A,b=sell->B; 813d4002b98SHong Zhang PetscErrorCode ierr; 814d4002b98SHong Zhang PetscInt s1,s2,s3; 815d4002b98SHong Zhang 816d4002b98SHong Zhang PetscFunctionBegin; 817d4002b98SHong Zhang ierr = MatGetLocalSize(mat,&s2,&s3);CHKERRQ(ierr); 818d4002b98SHong Zhang if (rr) { 819d4002b98SHong Zhang ierr = VecGetLocalSize(rr,&s1);CHKERRQ(ierr); 820d4002b98SHong Zhang if (s1!=s3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"right vector non-conforming local size"); 821d4002b98SHong Zhang /* Overlap communication with computation. */ 822d4002b98SHong Zhang ierr = VecScatterBegin(sell->Mvctx,rr,sell->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 823d4002b98SHong Zhang } 824d4002b98SHong Zhang if (ll) { 825d4002b98SHong Zhang ierr = VecGetLocalSize(ll,&s1);CHKERRQ(ierr); 826d4002b98SHong Zhang if (s1!=s2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"left vector non-conforming local size"); 827d4002b98SHong Zhang ierr = (*b->ops->diagonalscale)(b,ll,0);CHKERRQ(ierr); 828d4002b98SHong Zhang } 829d4002b98SHong Zhang /* scale the diagonal block */ 830d4002b98SHong Zhang ierr = (*a->ops->diagonalscale)(a,ll,rr);CHKERRQ(ierr); 831d4002b98SHong Zhang 832d4002b98SHong Zhang if (rr) { 833d4002b98SHong Zhang /* Do a scatter end and then right scale the off-diagonal block */ 834d4002b98SHong Zhang ierr = VecScatterEnd(sell->Mvctx,rr,sell->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 835d4002b98SHong Zhang ierr = (*b->ops->diagonalscale)(b,0,sell->lvec);CHKERRQ(ierr); 836d4002b98SHong Zhang } 837d4002b98SHong Zhang PetscFunctionReturn(0); 838d4002b98SHong Zhang } 839d4002b98SHong Zhang 840d4002b98SHong Zhang PetscErrorCode MatSetUnfactored_MPISELL(Mat A) 841d4002b98SHong Zhang { 842d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 843d4002b98SHong Zhang PetscErrorCode ierr; 844d4002b98SHong Zhang 845d4002b98SHong Zhang PetscFunctionBegin; 846d4002b98SHong Zhang ierr = MatSetUnfactored(a->A);CHKERRQ(ierr); 847d4002b98SHong Zhang PetscFunctionReturn(0); 848d4002b98SHong Zhang } 849d4002b98SHong Zhang 850d4002b98SHong Zhang PetscErrorCode MatEqual_MPISELL(Mat A,Mat B,PetscBool *flag) 851d4002b98SHong Zhang { 852d4002b98SHong Zhang Mat_MPISELL *matB=(Mat_MPISELL*)B->data,*matA=(Mat_MPISELL*)A->data; 853d4002b98SHong Zhang Mat a,b,c,d; 854d4002b98SHong Zhang PetscBool flg; 855d4002b98SHong Zhang PetscErrorCode ierr; 856d4002b98SHong Zhang 857d4002b98SHong Zhang PetscFunctionBegin; 858d4002b98SHong Zhang a = matA->A; b = matA->B; 859d4002b98SHong Zhang c = matB->A; d = matB->B; 860d4002b98SHong Zhang 861d4002b98SHong Zhang ierr = MatEqual(a,c,&flg);CHKERRQ(ierr); 862d4002b98SHong Zhang if (flg) { 863d4002b98SHong Zhang ierr = MatEqual(b,d,&flg);CHKERRQ(ierr); 864d4002b98SHong Zhang } 865d4002b98SHong Zhang ierr = MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 866d4002b98SHong Zhang PetscFunctionReturn(0); 867d4002b98SHong Zhang } 868d4002b98SHong Zhang 869d4002b98SHong Zhang PetscErrorCode MatCopy_MPISELL(Mat A,Mat B,MatStructure str) 870d4002b98SHong Zhang { 871d4002b98SHong Zhang PetscErrorCode ierr; 872d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 873d4002b98SHong Zhang Mat_MPISELL *b=(Mat_MPISELL*)B->data; 874d4002b98SHong Zhang 875d4002b98SHong Zhang PetscFunctionBegin; 876d4002b98SHong Zhang /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 877d4002b98SHong Zhang if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) { 878d4002b98SHong Zhang /* because of the column compression in the off-processor part of the matrix a->B, 879d4002b98SHong Zhang the number of columns in a->B and b->B may be different, hence we cannot call 880d4002b98SHong Zhang the MatCopy() directly on the two parts. If need be, we can provide a more 881d4002b98SHong Zhang efficient copy than the MatCopy_Basic() by first uncompressing the a->B matrices 882d4002b98SHong Zhang then copying the submatrices */ 883d4002b98SHong Zhang ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 884d4002b98SHong Zhang } else { 885d4002b98SHong Zhang ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr); 886d4002b98SHong Zhang ierr = MatCopy(a->B,b->B,str);CHKERRQ(ierr); 887d4002b98SHong Zhang } 888d4002b98SHong Zhang PetscFunctionReturn(0); 889d4002b98SHong Zhang } 890d4002b98SHong Zhang 891d4002b98SHong Zhang PetscErrorCode MatSetUp_MPISELL(Mat A) 892d4002b98SHong Zhang { 893d4002b98SHong Zhang PetscErrorCode ierr; 894d4002b98SHong Zhang 895d4002b98SHong Zhang PetscFunctionBegin; 896d4002b98SHong Zhang ierr = MatMPISELLSetPreallocation(A,PETSC_DEFAULT,0,PETSC_DEFAULT,0);CHKERRQ(ierr); 897d4002b98SHong Zhang PetscFunctionReturn(0); 898d4002b98SHong Zhang } 899d4002b98SHong Zhang 900d4002b98SHong Zhang 901d4002b98SHong Zhang extern PetscErrorCode MatConjugate_SeqSELL(Mat); 902d4002b98SHong Zhang 903d4002b98SHong Zhang PetscErrorCode MatConjugate_MPISELL(Mat mat) 904d4002b98SHong Zhang { 905d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX) 906d4002b98SHong Zhang PetscErrorCode ierr; 907d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 908d4002b98SHong Zhang 909d4002b98SHong Zhang PetscFunctionBegin; 910d4002b98SHong Zhang ierr = MatConjugate_SeqSELL(sell->A);CHKERRQ(ierr); 911d4002b98SHong Zhang ierr = MatConjugate_SeqSELL(sell->B);CHKERRQ(ierr); 912d4002b98SHong Zhang #else 913d4002b98SHong Zhang PetscFunctionBegin; 914d4002b98SHong Zhang #endif 915d4002b98SHong Zhang PetscFunctionReturn(0); 916d4002b98SHong Zhang } 917d4002b98SHong Zhang 918d4002b98SHong Zhang PetscErrorCode MatRealPart_MPISELL(Mat A) 919d4002b98SHong Zhang { 920d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 921d4002b98SHong Zhang PetscErrorCode ierr; 922d4002b98SHong Zhang 923d4002b98SHong Zhang PetscFunctionBegin; 924d4002b98SHong Zhang ierr = MatRealPart(a->A);CHKERRQ(ierr); 925d4002b98SHong Zhang ierr = MatRealPart(a->B);CHKERRQ(ierr); 926d4002b98SHong Zhang PetscFunctionReturn(0); 927d4002b98SHong Zhang } 928d4002b98SHong Zhang 929d4002b98SHong Zhang PetscErrorCode MatImaginaryPart_MPISELL(Mat A) 930d4002b98SHong Zhang { 931d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 932d4002b98SHong Zhang PetscErrorCode ierr; 933d4002b98SHong Zhang 934d4002b98SHong Zhang PetscFunctionBegin; 935d4002b98SHong Zhang ierr = MatImaginaryPart(a->A);CHKERRQ(ierr); 936d4002b98SHong Zhang ierr = MatImaginaryPart(a->B);CHKERRQ(ierr); 937d4002b98SHong Zhang PetscFunctionReturn(0); 938d4002b98SHong Zhang } 939d4002b98SHong Zhang 940d4002b98SHong Zhang PetscErrorCode MatInvertBlockDiagonal_MPISELL(Mat A,const PetscScalar **values) 941d4002b98SHong Zhang { 942d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 943d4002b98SHong Zhang PetscErrorCode ierr; 944d4002b98SHong Zhang 945d4002b98SHong Zhang PetscFunctionBegin; 946d4002b98SHong Zhang ierr = MatInvertBlockDiagonal(a->A,values);CHKERRQ(ierr); 947d4002b98SHong Zhang A->factorerrortype = a->A->factorerrortype; 948d4002b98SHong Zhang PetscFunctionReturn(0); 949d4002b98SHong Zhang } 950d4002b98SHong Zhang 951d4002b98SHong Zhang static PetscErrorCode MatSetRandom_MPISELL(Mat x,PetscRandom rctx) 952d4002b98SHong Zhang { 953d4002b98SHong Zhang PetscErrorCode ierr; 954d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)x->data; 955d4002b98SHong Zhang 956d4002b98SHong Zhang PetscFunctionBegin; 957d4002b98SHong Zhang ierr = MatSetRandom(sell->A,rctx);CHKERRQ(ierr); 958d4002b98SHong Zhang ierr = MatSetRandom(sell->B,rctx);CHKERRQ(ierr); 959d4002b98SHong Zhang ierr = MatAssemblyBegin(x,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 960d4002b98SHong Zhang ierr = MatAssemblyEnd(x,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 961d4002b98SHong Zhang PetscFunctionReturn(0); 962d4002b98SHong Zhang } 963d4002b98SHong Zhang 964d4002b98SHong Zhang PetscErrorCode MatSetFromOptions_MPISELL(PetscOptionItems *PetscOptionsObject,Mat A) 965d4002b98SHong Zhang { 966d4002b98SHong Zhang PetscErrorCode ierr; 967d4002b98SHong Zhang 968d4002b98SHong Zhang PetscFunctionBegin; 969d4002b98SHong Zhang ierr = PetscOptionsHead(PetscOptionsObject,"MPISELL options");CHKERRQ(ierr); 9700af67c1bSStefano Zampini ierr = PetscOptionsTail();CHKERRQ(ierr); 971d4002b98SHong Zhang PetscFunctionReturn(0); 972d4002b98SHong Zhang } 973d4002b98SHong Zhang 974d4002b98SHong Zhang PetscErrorCode MatShift_MPISELL(Mat Y,PetscScalar a) 975d4002b98SHong Zhang { 976d4002b98SHong Zhang PetscErrorCode ierr; 977d4002b98SHong Zhang Mat_MPISELL *msell=(Mat_MPISELL*)Y->data; 978d4002b98SHong Zhang Mat_SeqSELL *sell=(Mat_SeqSELL*)msell->A->data; 979d4002b98SHong Zhang 980d4002b98SHong Zhang PetscFunctionBegin; 981d4002b98SHong Zhang if (!Y->preallocated) { 982d4002b98SHong Zhang ierr = MatMPISELLSetPreallocation(Y,1,NULL,0,NULL);CHKERRQ(ierr); 983d4002b98SHong Zhang } else if (!sell->nz) { 984d4002b98SHong Zhang PetscInt nonew = sell->nonew; 985d4002b98SHong Zhang ierr = MatSeqSELLSetPreallocation(msell->A,1,NULL);CHKERRQ(ierr); 986d4002b98SHong Zhang sell->nonew = nonew; 987d4002b98SHong Zhang } 988d4002b98SHong Zhang ierr = MatShift_Basic(Y,a);CHKERRQ(ierr); 989d4002b98SHong Zhang PetscFunctionReturn(0); 990d4002b98SHong Zhang } 991d4002b98SHong Zhang 992d4002b98SHong Zhang PetscErrorCode MatMissingDiagonal_MPISELL(Mat A,PetscBool *missing,PetscInt *d) 993d4002b98SHong Zhang { 994d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 995d4002b98SHong Zhang PetscErrorCode ierr; 996d4002b98SHong Zhang 997d4002b98SHong Zhang PetscFunctionBegin; 998d4002b98SHong Zhang if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only works for square matrices"); 999d4002b98SHong Zhang ierr = MatMissingDiagonal(a->A,missing,d);CHKERRQ(ierr); 1000d4002b98SHong Zhang if (d) { 1001d4002b98SHong Zhang PetscInt rstart; 1002d4002b98SHong Zhang ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 1003d4002b98SHong Zhang *d += rstart; 1004d4002b98SHong Zhang 1005d4002b98SHong Zhang } 1006d4002b98SHong Zhang PetscFunctionReturn(0); 1007d4002b98SHong Zhang } 1008d4002b98SHong Zhang 1009d4002b98SHong Zhang PetscErrorCode MatGetDiagonalBlock_MPISELL(Mat A,Mat *a) 1010d4002b98SHong Zhang { 1011d4002b98SHong Zhang PetscFunctionBegin; 1012d4002b98SHong Zhang *a = ((Mat_MPISELL*)A->data)->A; 1013d4002b98SHong Zhang PetscFunctionReturn(0); 1014d4002b98SHong Zhang } 1015d4002b98SHong Zhang 1016d4002b98SHong Zhang /* -------------------------------------------------------------------*/ 1017d4002b98SHong Zhang static struct _MatOps MatOps_Values = {MatSetValues_MPISELL, 1018d4002b98SHong Zhang 0, 1019d4002b98SHong Zhang 0, 1020d4002b98SHong Zhang MatMult_MPISELL, 1021d4002b98SHong Zhang /* 4*/ MatMultAdd_MPISELL, 1022d4002b98SHong Zhang MatMultTranspose_MPISELL, 1023d4002b98SHong Zhang MatMultTransposeAdd_MPISELL, 1024d4002b98SHong Zhang 0, 1025d4002b98SHong Zhang 0, 1026d4002b98SHong Zhang 0, 1027d4002b98SHong Zhang /*10*/ 0, 1028d4002b98SHong Zhang 0, 1029d4002b98SHong Zhang 0, 1030d4002b98SHong Zhang MatSOR_MPISELL, 1031d4002b98SHong Zhang 0, 1032d4002b98SHong Zhang /*15*/ MatGetInfo_MPISELL, 1033d4002b98SHong Zhang MatEqual_MPISELL, 1034d4002b98SHong Zhang MatGetDiagonal_MPISELL, 1035d4002b98SHong Zhang MatDiagonalScale_MPISELL, 1036d4002b98SHong Zhang 0, 1037d4002b98SHong Zhang /*20*/ MatAssemblyBegin_MPISELL, 1038d4002b98SHong Zhang MatAssemblyEnd_MPISELL, 1039d4002b98SHong Zhang MatSetOption_MPISELL, 1040d4002b98SHong Zhang MatZeroEntries_MPISELL, 1041d4002b98SHong Zhang /*24*/ 0, 1042d4002b98SHong Zhang 0, 1043d4002b98SHong Zhang 0, 1044d4002b98SHong Zhang 0, 1045d4002b98SHong Zhang 0, 1046d4002b98SHong Zhang /*29*/ MatSetUp_MPISELL, 1047d4002b98SHong Zhang 0, 1048d4002b98SHong Zhang 0, 1049d4002b98SHong Zhang MatGetDiagonalBlock_MPISELL, 1050d4002b98SHong Zhang 0, 1051d4002b98SHong Zhang /*34*/ MatDuplicate_MPISELL, 1052d4002b98SHong Zhang 0, 1053d4002b98SHong Zhang 0, 1054d4002b98SHong Zhang 0, 1055d4002b98SHong Zhang 0, 1056d4002b98SHong Zhang /*39*/ 0, 1057d4002b98SHong Zhang 0, 1058d4002b98SHong Zhang 0, 1059d4002b98SHong Zhang MatGetValues_MPISELL, 1060d4002b98SHong Zhang MatCopy_MPISELL, 1061d4002b98SHong Zhang /*44*/ 0, 1062d4002b98SHong Zhang MatScale_MPISELL, 1063d4002b98SHong Zhang MatShift_MPISELL, 1064d4002b98SHong Zhang MatDiagonalSet_MPISELL, 1065d4002b98SHong Zhang 0, 1066d4002b98SHong Zhang /*49*/ MatSetRandom_MPISELL, 1067d4002b98SHong Zhang 0, 1068d4002b98SHong Zhang 0, 1069d4002b98SHong Zhang 0, 1070d4002b98SHong Zhang 0, 1071d4002b98SHong Zhang /*54*/ MatFDColoringCreate_MPIXAIJ, 1072d4002b98SHong Zhang 0, 1073d4002b98SHong Zhang MatSetUnfactored_MPISELL, 1074d4002b98SHong Zhang 0, 1075d4002b98SHong Zhang 0, 1076d4002b98SHong Zhang /*59*/ 0, 1077d4002b98SHong Zhang MatDestroy_MPISELL, 1078d4002b98SHong Zhang MatView_MPISELL, 1079d4002b98SHong Zhang 0, 1080d4002b98SHong Zhang 0, 1081d4002b98SHong Zhang /*64*/ 0, 1082d4002b98SHong Zhang 0, 1083d4002b98SHong Zhang 0, 1084d4002b98SHong Zhang 0, 1085d4002b98SHong Zhang 0, 1086d4002b98SHong Zhang /*69*/ 0, 1087d4002b98SHong Zhang 0, 1088d4002b98SHong Zhang 0, 1089d4002b98SHong Zhang 0, 1090d4002b98SHong Zhang 0, 1091d4002b98SHong Zhang 0, 1092d4002b98SHong Zhang /*75*/ MatFDColoringApply_AIJ, /* reuse AIJ function */ 1093d4002b98SHong Zhang MatSetFromOptions_MPISELL, 1094d4002b98SHong Zhang 0, 1095d4002b98SHong Zhang 0, 1096d4002b98SHong Zhang 0, 1097d4002b98SHong Zhang /*80*/ 0, 1098d4002b98SHong Zhang 0, 1099d4002b98SHong Zhang 0, 1100d4002b98SHong Zhang /*83*/ 0, 1101d4002b98SHong Zhang 0, 1102d4002b98SHong Zhang 0, 1103d4002b98SHong Zhang 0, 1104d4002b98SHong Zhang 0, 1105d4002b98SHong Zhang 0, 1106d4002b98SHong Zhang /*89*/ 0, 1107d4002b98SHong Zhang 0, 1108d4002b98SHong Zhang 0, 1109d4002b98SHong Zhang 0, 1110d4002b98SHong Zhang 0, 1111d4002b98SHong Zhang /*94*/ 0, 1112d4002b98SHong Zhang 0, 1113d4002b98SHong Zhang 0, 1114d4002b98SHong Zhang 0, 1115d4002b98SHong Zhang 0, 1116d4002b98SHong Zhang /*99*/ 0, 1117d4002b98SHong Zhang 0, 1118d4002b98SHong Zhang 0, 1119d4002b98SHong Zhang MatConjugate_MPISELL, 1120d4002b98SHong Zhang 0, 1121d4002b98SHong Zhang /*104*/0, 1122d4002b98SHong Zhang MatRealPart_MPISELL, 1123d4002b98SHong Zhang MatImaginaryPart_MPISELL, 1124d4002b98SHong Zhang 0, 1125d4002b98SHong Zhang 0, 1126d4002b98SHong Zhang /*109*/0, 1127d4002b98SHong Zhang 0, 1128d4002b98SHong Zhang 0, 1129d4002b98SHong Zhang 0, 1130d4002b98SHong Zhang MatMissingDiagonal_MPISELL, 1131d4002b98SHong Zhang /*114*/0, 1132d4002b98SHong Zhang 0, 1133d4002b98SHong Zhang MatGetGhosts_MPISELL, 1134d4002b98SHong Zhang 0, 1135d4002b98SHong Zhang 0, 1136d4002b98SHong Zhang /*119*/0, 1137d4002b98SHong Zhang 0, 1138d4002b98SHong Zhang 0, 1139d4002b98SHong Zhang 0, 1140d4002b98SHong Zhang 0, 1141d4002b98SHong Zhang /*124*/0, 1142d4002b98SHong Zhang 0, 1143d4002b98SHong Zhang MatInvertBlockDiagonal_MPISELL, 1144d4002b98SHong Zhang 0, 1145d4002b98SHong Zhang 0, 1146d4002b98SHong Zhang /*129*/0, 1147d4002b98SHong Zhang 0, 1148d4002b98SHong Zhang 0, 1149d4002b98SHong Zhang 0, 1150d4002b98SHong Zhang 0, 1151d4002b98SHong Zhang /*134*/0, 1152d4002b98SHong Zhang 0, 1153d4002b98SHong Zhang 0, 1154d4002b98SHong Zhang 0, 1155d4002b98SHong Zhang 0, 1156d4002b98SHong Zhang /*139*/0, 1157d4002b98SHong Zhang 0, 1158d4002b98SHong Zhang 0, 1159d4002b98SHong Zhang MatFDColoringSetUp_MPIXAIJ, 1160d4002b98SHong Zhang 0, 1161d4002b98SHong Zhang /*144*/0 1162d4002b98SHong Zhang }; 1163d4002b98SHong Zhang 1164d4002b98SHong Zhang /* ----------------------------------------------------------------------------------------*/ 1165d4002b98SHong Zhang 1166d4002b98SHong Zhang PetscErrorCode MatStoreValues_MPISELL(Mat mat) 1167d4002b98SHong Zhang { 1168d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 1169d4002b98SHong Zhang PetscErrorCode ierr; 1170d4002b98SHong Zhang 1171d4002b98SHong Zhang PetscFunctionBegin; 1172d4002b98SHong Zhang ierr = MatStoreValues(sell->A);CHKERRQ(ierr); 1173d4002b98SHong Zhang ierr = MatStoreValues(sell->B);CHKERRQ(ierr); 1174d4002b98SHong Zhang PetscFunctionReturn(0); 1175d4002b98SHong Zhang } 1176d4002b98SHong Zhang 1177d4002b98SHong Zhang PetscErrorCode MatRetrieveValues_MPISELL(Mat mat) 1178d4002b98SHong Zhang { 1179d4002b98SHong Zhang Mat_MPISELL *sell=(Mat_MPISELL*)mat->data; 1180d4002b98SHong Zhang PetscErrorCode ierr; 1181d4002b98SHong Zhang 1182d4002b98SHong Zhang PetscFunctionBegin; 1183d4002b98SHong Zhang ierr = MatRetrieveValues(sell->A);CHKERRQ(ierr); 1184d4002b98SHong Zhang ierr = MatRetrieveValues(sell->B);CHKERRQ(ierr); 1185d4002b98SHong Zhang PetscFunctionReturn(0); 1186d4002b98SHong Zhang } 1187d4002b98SHong Zhang 1188d4002b98SHong Zhang PetscErrorCode MatMPISELLSetPreallocation_MPISELL(Mat B,PetscInt d_rlenmax,const PetscInt d_rlen[],PetscInt o_rlenmax,const PetscInt o_rlen[]) 1189d4002b98SHong Zhang { 1190d4002b98SHong Zhang Mat_MPISELL *b; 1191d4002b98SHong Zhang PetscErrorCode ierr; 1192d4002b98SHong Zhang 1193d4002b98SHong Zhang PetscFunctionBegin; 1194d4002b98SHong Zhang ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 1195d4002b98SHong Zhang ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 1196d4002b98SHong Zhang b = (Mat_MPISELL*)B->data; 1197d4002b98SHong Zhang 1198d4002b98SHong Zhang if (!B->preallocated) { 1199d4002b98SHong Zhang /* Explicitly create 2 MATSEQSELL matrices. */ 1200d4002b98SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&b->A);CHKERRQ(ierr); 1201d4002b98SHong Zhang ierr = MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);CHKERRQ(ierr); 1202d4002b98SHong Zhang ierr = MatSetBlockSizesFromMats(b->A,B,B);CHKERRQ(ierr); 1203d4002b98SHong Zhang ierr = MatSetType(b->A,MATSEQSELL);CHKERRQ(ierr); 1204d4002b98SHong Zhang ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)b->A);CHKERRQ(ierr); 1205d4002b98SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&b->B);CHKERRQ(ierr); 1206d4002b98SHong Zhang ierr = MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N);CHKERRQ(ierr); 1207d4002b98SHong Zhang ierr = MatSetBlockSizesFromMats(b->B,B,B);CHKERRQ(ierr); 1208d4002b98SHong Zhang ierr = MatSetType(b->B,MATSEQSELL);CHKERRQ(ierr); 1209d4002b98SHong Zhang ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)b->B);CHKERRQ(ierr); 1210d4002b98SHong Zhang } 1211d4002b98SHong Zhang 1212d4002b98SHong Zhang ierr = MatSeqSELLSetPreallocation(b->A,d_rlenmax,d_rlen);CHKERRQ(ierr); 1213d4002b98SHong Zhang ierr = MatSeqSELLSetPreallocation(b->B,o_rlenmax,o_rlen);CHKERRQ(ierr); 1214d4002b98SHong Zhang B->preallocated = PETSC_TRUE; 1215d4002b98SHong Zhang B->was_assembled = PETSC_FALSE; 1216d4002b98SHong Zhang /* 1217d4002b98SHong Zhang critical for MatAssemblyEnd to work. 1218d4002b98SHong Zhang MatAssemblyBegin checks it to set up was_assembled 1219d4002b98SHong Zhang and MatAssemblyEnd checks was_assembled to determine whether to build garray 1220d4002b98SHong Zhang */ 1221d4002b98SHong Zhang B->assembled = PETSC_FALSE; 1222d4002b98SHong Zhang PetscFunctionReturn(0); 1223d4002b98SHong Zhang } 1224d4002b98SHong Zhang 1225d4002b98SHong Zhang PetscErrorCode MatDuplicate_MPISELL(Mat matin,MatDuplicateOption cpvalues,Mat *newmat) 1226d4002b98SHong Zhang { 1227d4002b98SHong Zhang Mat mat; 1228d4002b98SHong Zhang Mat_MPISELL *a,*oldmat=(Mat_MPISELL*)matin->data; 1229d4002b98SHong Zhang PetscErrorCode ierr; 1230d4002b98SHong Zhang 1231d4002b98SHong Zhang PetscFunctionBegin; 1232d4002b98SHong Zhang *newmat = 0; 1233d4002b98SHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)matin),&mat);CHKERRQ(ierr); 1234d4002b98SHong Zhang ierr = MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N);CHKERRQ(ierr); 1235d4002b98SHong Zhang ierr = MatSetBlockSizesFromMats(mat,matin,matin);CHKERRQ(ierr); 1236d4002b98SHong Zhang ierr = MatSetType(mat,((PetscObject)matin)->type_name);CHKERRQ(ierr); 1237d4002b98SHong Zhang a = (Mat_MPISELL*)mat->data; 1238d4002b98SHong Zhang 1239d4002b98SHong Zhang mat->factortype = matin->factortype; 1240d4002b98SHong Zhang mat->assembled = PETSC_TRUE; 1241d4002b98SHong Zhang mat->insertmode = NOT_SET_VALUES; 1242d4002b98SHong Zhang mat->preallocated = PETSC_TRUE; 1243d4002b98SHong Zhang 1244d4002b98SHong Zhang a->size = oldmat->size; 1245d4002b98SHong Zhang a->rank = oldmat->rank; 1246d4002b98SHong Zhang a->donotstash = oldmat->donotstash; 1247d4002b98SHong Zhang a->roworiented = oldmat->roworiented; 1248d4002b98SHong Zhang a->rowindices = 0; 1249d4002b98SHong Zhang a->rowvalues = 0; 1250d4002b98SHong Zhang a->getrowactive = PETSC_FALSE; 1251d4002b98SHong Zhang 1252d4002b98SHong Zhang ierr = PetscLayoutReference(matin->rmap,&mat->rmap);CHKERRQ(ierr); 1253d4002b98SHong Zhang ierr = PetscLayoutReference(matin->cmap,&mat->cmap);CHKERRQ(ierr); 1254d4002b98SHong Zhang 1255d4002b98SHong Zhang if (oldmat->colmap) { 1256d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE) 1257d4002b98SHong Zhang ierr = PetscTableCreateCopy(oldmat->colmap,&a->colmap);CHKERRQ(ierr); 1258d4002b98SHong Zhang #else 1259d4002b98SHong Zhang ierr = PetscMalloc1(mat->cmap->N,&a->colmap);CHKERRQ(ierr); 1260d4002b98SHong Zhang ierr = PetscLogObjectMemory((PetscObject)mat,(mat->cmap->N)*sizeof(PetscInt));CHKERRQ(ierr); 1261d4002b98SHong Zhang ierr = PetscMemcpy(a->colmap,oldmat->colmap,(mat->cmap->N)*sizeof(PetscInt));CHKERRQ(ierr); 1262d4002b98SHong Zhang #endif 1263d4002b98SHong Zhang } else a->colmap = 0; 1264d4002b98SHong Zhang if (oldmat->garray) { 1265d4002b98SHong Zhang PetscInt len; 1266d4002b98SHong Zhang len = oldmat->B->cmap->n; 1267d4002b98SHong Zhang ierr = PetscMalloc1(len+1,&a->garray);CHKERRQ(ierr); 1268d4002b98SHong Zhang ierr = PetscLogObjectMemory((PetscObject)mat,len*sizeof(PetscInt));CHKERRQ(ierr); 1269d4002b98SHong Zhang if (len) { ierr = PetscMemcpy(a->garray,oldmat->garray,len*sizeof(PetscInt));CHKERRQ(ierr); } 1270d4002b98SHong Zhang } else a->garray = 0; 1271d4002b98SHong Zhang 1272d4002b98SHong Zhang ierr = VecDuplicate(oldmat->lvec,&a->lvec);CHKERRQ(ierr); 1273d4002b98SHong Zhang ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->lvec);CHKERRQ(ierr); 1274d4002b98SHong Zhang ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx);CHKERRQ(ierr); 1275d4002b98SHong Zhang ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->Mvctx);CHKERRQ(ierr); 1276d4002b98SHong Zhang ierr = MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr); 1277d4002b98SHong Zhang ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);CHKERRQ(ierr); 1278d4002b98SHong Zhang ierr = MatDuplicate(oldmat->B,cpvalues,&a->B);CHKERRQ(ierr); 1279d4002b98SHong Zhang ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)a->B);CHKERRQ(ierr); 1280d4002b98SHong Zhang ierr = PetscFunctionListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist);CHKERRQ(ierr); 1281d4002b98SHong Zhang *newmat = mat; 1282d4002b98SHong Zhang PetscFunctionReturn(0); 1283d4002b98SHong Zhang } 1284d4002b98SHong Zhang 1285d4002b98SHong Zhang /*@C 1286d4002b98SHong Zhang MatMPISELLSetPreallocation - Preallocates memory for a sparse parallel matrix in sell format. 1287d4002b98SHong Zhang For good matrix assembly performance the user should preallocate the matrix storage by 1288d4002b98SHong Zhang setting the parameters d_nz (or d_nnz) and o_nz (or o_nnz). 1289d4002b98SHong Zhang 1290d4002b98SHong Zhang Collective on MPI_Comm 1291d4002b98SHong Zhang 1292d4002b98SHong Zhang Input Parameters: 1293d4002b98SHong Zhang + B - the matrix 1294d4002b98SHong Zhang . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 1295d4002b98SHong Zhang (same value is used for all local rows) 1296d4002b98SHong Zhang . d_nnz - array containing the number of nonzeros in the various rows of the 1297d4002b98SHong Zhang DIAGONAL portion of the local submatrix (possibly different for each row) 1298d4002b98SHong Zhang or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure. 1299d4002b98SHong Zhang The size of this array is equal to the number of local rows, i.e 'm'. 1300d4002b98SHong Zhang For matrices that will be factored, you must leave room for (and set) 1301d4002b98SHong Zhang the diagonal entry even if it is zero. 1302d4002b98SHong Zhang . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 1303d4002b98SHong Zhang submatrix (same value is used for all local rows). 1304d4002b98SHong Zhang - o_nnz - array containing the number of nonzeros in the various rows of the 1305d4002b98SHong Zhang OFF-DIAGONAL portion of the local submatrix (possibly different for 1306d4002b98SHong Zhang each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero 1307d4002b98SHong Zhang structure. The size of this array is equal to the number 1308d4002b98SHong Zhang of local rows, i.e 'm'. 1309d4002b98SHong Zhang 1310d4002b98SHong Zhang If the *_nnz parameter is given then the *_nz parameter is ignored 1311d4002b98SHong Zhang 1312d4002b98SHong Zhang The stored row and column indices begin with zero. 1313d4002b98SHong Zhang 1314d4002b98SHong Zhang The parallel matrix is partitioned such that the first m0 rows belong to 1315d4002b98SHong Zhang process 0, the next m1 rows belong to process 1, the next m2 rows belong 1316d4002b98SHong Zhang to process 2 etc.. where m0,m1,m2... are the input parameter 'm'. 1317d4002b98SHong Zhang 1318d4002b98SHong Zhang The DIAGONAL portion of the local submatrix of a processor can be defined 1319d4002b98SHong Zhang as the submatrix which is obtained by extraction the part corresponding to 1320d4002b98SHong Zhang the rows r1-r2 and columns c1-c2 of the global matrix, where r1 is the 1321d4002b98SHong Zhang first row that belongs to the processor, r2 is the last row belonging to 1322d4002b98SHong Zhang the this processor, and c1-c2 is range of indices of the local part of a 1323d4002b98SHong Zhang vector suitable for applying the matrix to. This is an mxn matrix. In the 1324d4002b98SHong Zhang common case of a square matrix, the row and column ranges are the same and 1325d4002b98SHong Zhang the DIAGONAL part is also square. The remaining portion of the local 1326d4002b98SHong Zhang submatrix (mxN) constitute the OFF-DIAGONAL portion. 1327d4002b98SHong Zhang 1328d4002b98SHong Zhang If o_nnz, d_nnz are specified, then o_nz, and d_nz are ignored. 1329d4002b98SHong Zhang 1330d4002b98SHong Zhang You can call MatGetInfo() to get information on how effective the preallocation was; 1331d4002b98SHong Zhang for example the fields mallocs,nz_allocated,nz_used,nz_unneeded; 1332d4002b98SHong Zhang You can also run with the option -info and look for messages with the string 1333d4002b98SHong Zhang malloc in them to see if additional memory allocation was needed. 1334d4002b98SHong Zhang 1335d4002b98SHong Zhang Example usage: 1336d4002b98SHong Zhang 1337d4002b98SHong Zhang Consider the following 8x8 matrix with 34 non-zero values, that is 1338d4002b98SHong Zhang assembled across 3 processors. Lets assume that proc0 owns 3 rows, 1339d4002b98SHong Zhang proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown 1340d4002b98SHong Zhang as follows: 1341d4002b98SHong Zhang 1342d4002b98SHong Zhang .vb 1343d4002b98SHong Zhang 1 2 0 | 0 3 0 | 0 4 1344d4002b98SHong Zhang Proc0 0 5 6 | 7 0 0 | 8 0 1345d4002b98SHong Zhang 9 0 10 | 11 0 0 | 12 0 1346d4002b98SHong Zhang ------------------------------------- 1347d4002b98SHong Zhang 13 0 14 | 15 16 17 | 0 0 1348d4002b98SHong Zhang Proc1 0 18 0 | 19 20 21 | 0 0 1349d4002b98SHong Zhang 0 0 0 | 22 23 0 | 24 0 1350d4002b98SHong Zhang ------------------------------------- 1351d4002b98SHong Zhang Proc2 25 26 27 | 0 0 28 | 29 0 1352d4002b98SHong Zhang 30 0 0 | 31 32 33 | 0 34 1353d4002b98SHong Zhang .ve 1354d4002b98SHong Zhang 1355d4002b98SHong Zhang This can be represented as a collection of submatrices as: 1356d4002b98SHong Zhang 1357d4002b98SHong Zhang .vb 1358d4002b98SHong Zhang A B C 1359d4002b98SHong Zhang D E F 1360d4002b98SHong Zhang G H I 1361d4002b98SHong Zhang .ve 1362d4002b98SHong Zhang 1363d4002b98SHong Zhang Where the submatrices A,B,C are owned by proc0, D,E,F are 1364d4002b98SHong Zhang owned by proc1, G,H,I are owned by proc2. 1365d4002b98SHong Zhang 1366d4002b98SHong Zhang The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively. 1367d4002b98SHong Zhang The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively. 1368d4002b98SHong Zhang The 'M','N' parameters are 8,8, and have the same values on all procs. 1369d4002b98SHong Zhang 1370d4002b98SHong Zhang The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are 1371d4002b98SHong Zhang submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices 1372d4002b98SHong Zhang corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively. 1373d4002b98SHong Zhang Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL 1374d4002b98SHong Zhang part as SeqSELL matrices. for eg: proc1 will store [E] as a SeqSELL 1375d4002b98SHong Zhang matrix, ans [DF] as another SeqSELL matrix. 1376d4002b98SHong Zhang 1377d4002b98SHong Zhang When d_nz, o_nz parameters are specified, d_nz storage elements are 1378d4002b98SHong Zhang allocated for every row of the local diagonal submatrix, and o_nz 1379d4002b98SHong Zhang storage locations are allocated for every row of the OFF-DIAGONAL submat. 1380d4002b98SHong Zhang One way to choose d_nz and o_nz is to use the max nonzerors per local 1381d4002b98SHong Zhang rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices. 1382d4002b98SHong Zhang In this case, the values of d_nz,o_nz are: 1383d4002b98SHong Zhang .vb 1384d4002b98SHong Zhang proc0 : dnz = 2, o_nz = 2 1385d4002b98SHong Zhang proc1 : dnz = 3, o_nz = 2 1386d4002b98SHong Zhang proc2 : dnz = 1, o_nz = 4 1387d4002b98SHong Zhang .ve 1388d4002b98SHong Zhang We are allocating m*(d_nz+o_nz) storage locations for every proc. This 1389d4002b98SHong Zhang translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10 1390d4002b98SHong Zhang for proc3. i.e we are using 12+15+10=37 storage locations to store 1391d4002b98SHong Zhang 34 values. 1392d4002b98SHong Zhang 1393d4002b98SHong Zhang When d_nnz, o_nnz parameters are specified, the storage is specified 1394d4002b98SHong Zhang for every row, coresponding to both DIAGONAL and OFF-DIAGONAL submatrices. 1395d4002b98SHong Zhang In the above case the values for d_nnz,o_nnz are: 1396d4002b98SHong Zhang .vb 1397d4002b98SHong Zhang proc0: d_nnz = [2,2,2] and o_nnz = [2,2,2] 1398d4002b98SHong Zhang proc1: d_nnz = [3,3,2] and o_nnz = [2,1,1] 1399d4002b98SHong Zhang proc2: d_nnz = [1,1] and o_nnz = [4,4] 1400d4002b98SHong Zhang .ve 1401d4002b98SHong Zhang Here the space allocated is according to nz (or maximum values in the nnz 1402d4002b98SHong Zhang if nnz is provided) for DIAGONAL and OFF-DIAGONAL submatrices, i.e (2+2+3+2)*3+(1+4)*2=37 1403d4002b98SHong Zhang 1404d4002b98SHong Zhang Level: intermediate 1405d4002b98SHong Zhang 1406d4002b98SHong Zhang .keywords: matrix, sell, sparse, parallel 1407d4002b98SHong Zhang 1408d4002b98SHong Zhang .seealso: MatCreate(), MatCreateSeqSELL(), MatSetValues(), MatCreatesell(), 1409d4002b98SHong Zhang MATMPISELL, MatGetInfo(), PetscSplitOwnership() 1410d4002b98SHong Zhang @*/ 1411d4002b98SHong Zhang PetscErrorCode MatMPISELLSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 1412d4002b98SHong Zhang { 1413d4002b98SHong Zhang PetscErrorCode ierr; 1414d4002b98SHong Zhang 1415d4002b98SHong Zhang PetscFunctionBegin; 1416d4002b98SHong Zhang PetscValidHeaderSpecific(B,MAT_CLASSID,1); 1417d4002b98SHong Zhang PetscValidType(B,1); 1418d4002b98SHong Zhang ierr = PetscTryMethod(B,"MatMPISELLSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 1419d4002b98SHong Zhang PetscFunctionReturn(0); 1420d4002b98SHong Zhang } 1421d4002b98SHong Zhang 1422d4002b98SHong Zhang /*@C 1423d4002b98SHong Zhang MatCreateSELL - Creates a sparse parallel matrix in SELL format. 1424d4002b98SHong Zhang 1425d4002b98SHong Zhang Collective on MPI_Comm 1426d4002b98SHong Zhang 1427d4002b98SHong Zhang Input Parameters: 1428d4002b98SHong Zhang + comm - MPI communicator 1429d4002b98SHong Zhang . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 1430d4002b98SHong Zhang This value should be the same as the local size used in creating the 1431d4002b98SHong Zhang y vector for the matrix-vector product y = Ax. 1432d4002b98SHong Zhang . n - This value should be the same as the local size used in creating the 1433d4002b98SHong Zhang x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 1434d4002b98SHong Zhang calculated if N is given) For square matrices n is almost always m. 1435d4002b98SHong Zhang . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 1436d4002b98SHong Zhang . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 1437d4002b98SHong Zhang . d_rlenmax - max number of nonzeros per row in DIAGONAL portion of local submatrix 1438d4002b98SHong Zhang (same value is used for all local rows) 1439d4002b98SHong Zhang . d_rlen - array containing the number of nonzeros in the various rows of the 1440d4002b98SHong Zhang DIAGONAL portion of the local submatrix (possibly different for each row) 1441d4002b98SHong Zhang or NULL, if d_rlenmax is used to specify the nonzero structure. 1442d4002b98SHong Zhang The size of this array is equal to the number of local rows, i.e 'm'. 1443d4002b98SHong Zhang . o_rlenmax - max number of nonzeros per row in the OFF-DIAGONAL portion of local 1444d4002b98SHong Zhang submatrix (same value is used for all local rows). 1445d4002b98SHong Zhang - o_rlen - array containing the number of nonzeros in the various rows of the 1446d4002b98SHong Zhang OFF-DIAGONAL portion of the local submatrix (possibly different for 1447d4002b98SHong Zhang each row) or NULL, if o_rlenmax is used to specify the nonzero 1448d4002b98SHong Zhang structure. The size of this array is equal to the number 1449d4002b98SHong Zhang of local rows, i.e 'm'. 1450d4002b98SHong Zhang 1451d4002b98SHong Zhang Output Parameter: 1452d4002b98SHong Zhang . A - the matrix 1453d4002b98SHong Zhang 1454d4002b98SHong Zhang It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(), 1455d4002b98SHong Zhang MatXXXXSetPreallocation() paradgm instead of this routine directly. 1456d4002b98SHong Zhang [MatXXXXSetPreallocation() is, for example, MatSeqSELLSetPreallocation] 1457d4002b98SHong Zhang 1458d4002b98SHong Zhang Notes: 1459d4002b98SHong Zhang If the *_rlen parameter is given then the *_rlenmax parameter is ignored 1460d4002b98SHong Zhang 1461d4002b98SHong Zhang m,n,M,N parameters specify the size of the matrix, and its partitioning across 1462d4002b98SHong Zhang processors, while d_rlenmax,d_rlen,o_rlenmax,o_rlen parameters specify the approximate 1463d4002b98SHong Zhang storage requirements for this matrix. 1464d4002b98SHong Zhang 1465d4002b98SHong Zhang If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one 1466d4002b98SHong Zhang processor than it must be used on all processors that share the object for 1467d4002b98SHong Zhang that argument. 1468d4002b98SHong Zhang 1469d4002b98SHong Zhang The user MUST specify either the local or global matrix dimensions 1470d4002b98SHong Zhang (possibly both). 1471d4002b98SHong Zhang 1472d4002b98SHong Zhang The parallel matrix is partitioned across processors such that the 1473d4002b98SHong Zhang first m0 rows belong to process 0, the next m1 rows belong to 1474d4002b98SHong Zhang process 1, the next m2 rows belong to process 2 etc.. where 1475d4002b98SHong Zhang m0,m1,m2,.. are the input parameter 'm'. i.e each processor stores 1476d4002b98SHong Zhang values corresponding to [m x N] submatrix. 1477d4002b98SHong Zhang 1478d4002b98SHong Zhang The columns are logically partitioned with the n0 columns belonging 1479d4002b98SHong Zhang to 0th partition, the next n1 columns belonging to the next 1480d4002b98SHong Zhang partition etc.. where n0,n1,n2... are the input parameter 'n'. 1481d4002b98SHong Zhang 1482d4002b98SHong Zhang The DIAGONAL portion of the local submatrix on any given processor 1483d4002b98SHong Zhang is the submatrix corresponding to the rows and columns m,n 1484d4002b98SHong Zhang corresponding to the given processor. i.e diagonal matrix on 1485d4002b98SHong Zhang process 0 is [m0 x n0], diagonal matrix on process 1 is [m1 x n1] 1486d4002b98SHong Zhang etc. The remaining portion of the local submatrix [m x (N-n)] 1487d4002b98SHong Zhang constitute the OFF-DIAGONAL portion. The example below better 1488d4002b98SHong Zhang illustrates this concept. 1489d4002b98SHong Zhang 1490d4002b98SHong Zhang For a square global matrix we define each processor's diagonal portion 1491d4002b98SHong Zhang to be its local rows and the corresponding columns (a square submatrix); 1492d4002b98SHong Zhang each processor's off-diagonal portion encompasses the remainder of the 1493d4002b98SHong Zhang local matrix (a rectangular submatrix). 1494d4002b98SHong Zhang 1495d4002b98SHong Zhang If o_rlen, d_rlen are specified, then o_rlenmax, and d_rlenmax are ignored. 1496d4002b98SHong Zhang 1497d4002b98SHong Zhang When calling this routine with a single process communicator, a matrix of 1498d4002b98SHong Zhang type SEQSELL is returned. If a matrix of type MATMPISELL is desired for this 1499d4002b98SHong Zhang type of communicator, use the construction mechanism: 1500d4002b98SHong Zhang MatCreate(...,&A); MatSetType(A,MATMPISELL); MatSetSizes(A, m,n,M,N); MatMPISELLSetPreallocation(A,...); 1501d4002b98SHong Zhang 1502d4002b98SHong Zhang Options Database Keys: 1503d4002b98SHong Zhang - -mat_sell_oneindex - Internally use indexing starting at 1 1504d4002b98SHong Zhang rather than 0. Note that when calling MatSetValues(), 1505d4002b98SHong Zhang the user still MUST index entries starting at 0! 1506d4002b98SHong Zhang 1507d4002b98SHong Zhang 1508d4002b98SHong Zhang Example usage: 1509d4002b98SHong Zhang 1510d4002b98SHong Zhang Consider the following 8x8 matrix with 34 non-zero values, that is 1511d4002b98SHong Zhang assembled across 3 processors. Lets assume that proc0 owns 3 rows, 1512d4002b98SHong Zhang proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown 1513d4002b98SHong Zhang as follows: 1514d4002b98SHong Zhang 1515d4002b98SHong Zhang .vb 1516d4002b98SHong Zhang 1 2 0 | 0 3 0 | 0 4 1517d4002b98SHong Zhang Proc0 0 5 6 | 7 0 0 | 8 0 1518d4002b98SHong Zhang 9 0 10 | 11 0 0 | 12 0 1519d4002b98SHong Zhang ------------------------------------- 1520d4002b98SHong Zhang 13 0 14 | 15 16 17 | 0 0 1521d4002b98SHong Zhang Proc1 0 18 0 | 19 20 21 | 0 0 1522d4002b98SHong Zhang 0 0 0 | 22 23 0 | 24 0 1523d4002b98SHong Zhang ------------------------------------- 1524d4002b98SHong Zhang Proc2 25 26 27 | 0 0 28 | 29 0 1525d4002b98SHong Zhang 30 0 0 | 31 32 33 | 0 34 1526d4002b98SHong Zhang .ve 1527d4002b98SHong Zhang 1528d4002b98SHong Zhang This can be represented as a collection of submatrices as: 1529d4002b98SHong Zhang 1530d4002b98SHong Zhang .vb 1531d4002b98SHong Zhang A B C 1532d4002b98SHong Zhang D E F 1533d4002b98SHong Zhang G H I 1534d4002b98SHong Zhang .ve 1535d4002b98SHong Zhang 1536d4002b98SHong Zhang Where the submatrices A,B,C are owned by proc0, D,E,F are 1537d4002b98SHong Zhang owned by proc1, G,H,I are owned by proc2. 1538d4002b98SHong Zhang 1539d4002b98SHong Zhang The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively. 1540d4002b98SHong Zhang The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively. 1541d4002b98SHong Zhang The 'M','N' parameters are 8,8, and have the same values on all procs. 1542d4002b98SHong Zhang 1543d4002b98SHong Zhang The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are 1544d4002b98SHong Zhang submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices 1545d4002b98SHong Zhang corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively. 1546d4002b98SHong Zhang Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL 1547d4002b98SHong Zhang part as SeqSELL matrices. for eg: proc1 will store [E] as a SeqSELL 1548d4002b98SHong Zhang matrix, ans [DF] as another SeqSELL matrix. 1549d4002b98SHong Zhang 1550d4002b98SHong Zhang When d_rlenmax, o_rlenmax parameters are specified, d_rlenmax storage elements are 1551d4002b98SHong Zhang allocated for every row of the local diagonal submatrix, and o_rlenmax 1552d4002b98SHong Zhang storage locations are allocated for every row of the OFF-DIAGONAL submat. 1553d4002b98SHong Zhang One way to choose d_rlenmax and o_rlenmax is to use the max nonzerors per local 1554d4002b98SHong Zhang rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices. 1555d4002b98SHong Zhang In this case, the values of d_rlenmax,o_rlenmax are: 1556d4002b98SHong Zhang .vb 1557d4002b98SHong Zhang proc0 : d_rlenmax = 2, o_rlenmax = 2 1558d4002b98SHong Zhang proc1 : d_rlenmax = 3, o_rlenmax = 2 1559d4002b98SHong Zhang proc2 : d_rlenmax = 1, o_rlenmax = 4 1560d4002b98SHong Zhang .ve 1561d4002b98SHong Zhang We are allocating m*(d_rlenmax+o_rlenmax) storage locations for every proc. This 1562d4002b98SHong Zhang translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10 1563d4002b98SHong Zhang for proc3. i.e we are using 12+15+10=37 storage locations to store 1564d4002b98SHong Zhang 34 values. 1565d4002b98SHong Zhang 1566d4002b98SHong Zhang When d_rlen, o_rlen parameters are specified, the storage is specified 1567d4002b98SHong Zhang for every row, coresponding to both DIAGONAL and OFF-DIAGONAL submatrices. 1568d4002b98SHong Zhang In the above case the values for d_nnz,o_nnz are: 1569d4002b98SHong Zhang .vb 1570d4002b98SHong Zhang proc0: d_nnz = [2,2,2] and o_nnz = [2,2,2] 1571d4002b98SHong Zhang proc1: d_nnz = [3,3,2] and o_nnz = [2,1,1] 1572d4002b98SHong Zhang proc2: d_nnz = [1,1] and o_nnz = [4,4] 1573d4002b98SHong Zhang .ve 1574d4002b98SHong Zhang Here the space allocated is still 37 though there are 34 nonzeros because 1575d4002b98SHong Zhang the allocation is always done according to rlenmax. 1576d4002b98SHong Zhang 1577d4002b98SHong Zhang Level: intermediate 1578d4002b98SHong Zhang 1579d4002b98SHong Zhang .keywords: matrix, sell, sparse, parallel 1580d4002b98SHong Zhang 1581d4002b98SHong Zhang .seealso: MatCreate(), MatCreateSeqSELL(), MatSetValues(), MatMPISELLSetPreallocation(), MatMPISELLSetPreallocationSELL(), 1582d4002b98SHong Zhang MATMPISELL, MatCreateMPISELLWithArrays() 1583d4002b98SHong Zhang @*/ 1584d4002b98SHong Zhang PetscErrorCode MatCreateSELL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt d_rlenmax,const PetscInt d_rlen[],PetscInt o_rlenmax,const PetscInt o_rlen[],Mat *A) 1585d4002b98SHong Zhang { 1586d4002b98SHong Zhang PetscErrorCode ierr; 1587d4002b98SHong Zhang PetscMPIInt size; 1588d4002b98SHong Zhang 1589d4002b98SHong Zhang PetscFunctionBegin; 1590d4002b98SHong Zhang ierr = MatCreate(comm,A);CHKERRQ(ierr); 1591d4002b98SHong Zhang ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1592d4002b98SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1593d4002b98SHong Zhang if (size > 1) { 1594d4002b98SHong Zhang ierr = MatSetType(*A,MATMPISELL);CHKERRQ(ierr); 1595d4002b98SHong Zhang ierr = MatMPISELLSetPreallocation(*A,d_rlenmax,d_rlen,o_rlenmax,o_rlen);CHKERRQ(ierr); 1596d4002b98SHong Zhang } else { 1597d4002b98SHong Zhang ierr = MatSetType(*A,MATSEQSELL);CHKERRQ(ierr); 1598d4002b98SHong Zhang ierr = MatSeqSELLSetPreallocation(*A,d_rlenmax,d_rlen);CHKERRQ(ierr); 1599d4002b98SHong Zhang } 1600d4002b98SHong Zhang PetscFunctionReturn(0); 1601d4002b98SHong Zhang } 1602d4002b98SHong Zhang 1603d4002b98SHong Zhang PetscErrorCode MatMPISELLGetSeqSELL(Mat A,Mat *Ad,Mat *Ao,const PetscInt *colmap[]) 1604d4002b98SHong Zhang { 1605d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 1606d4002b98SHong Zhang PetscBool flg; 1607d4002b98SHong Zhang PetscErrorCode ierr; 1608d4002b98SHong Zhang 1609d4002b98SHong Zhang PetscFunctionBegin; 1610d4002b98SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)A,MATMPISELL,&flg);CHKERRQ(ierr); 1611d4002b98SHong Zhang if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"This function requires a MATMPISELL matrix as input"); 1612d4002b98SHong Zhang if (Ad) *Ad = a->A; 1613d4002b98SHong Zhang if (Ao) *Ao = a->B; 1614d4002b98SHong Zhang if (colmap) *colmap = a->garray; 1615d4002b98SHong Zhang PetscFunctionReturn(0); 1616d4002b98SHong Zhang } 1617d4002b98SHong Zhang 1618d4002b98SHong Zhang /*@C 1619d4002b98SHong Zhang MatMPISELLGetLocalMatCondensed - Creates a SeqSELL matrix from an MATMPISELL matrix by taking all its local rows and NON-ZERO columns 1620d4002b98SHong Zhang 1621d4002b98SHong Zhang Not Collective 1622d4002b98SHong Zhang 1623d4002b98SHong Zhang Input Parameters: 1624d4002b98SHong Zhang + A - the matrix 1625d4002b98SHong Zhang . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 1626d4002b98SHong Zhang - row, col - index sets of rows and columns to extract (or NULL) 1627d4002b98SHong Zhang 1628d4002b98SHong Zhang Output Parameter: 1629d4002b98SHong Zhang . A_loc - the local sequential matrix generated 1630d4002b98SHong Zhang 1631d4002b98SHong Zhang Level: developer 1632d4002b98SHong Zhang 1633d4002b98SHong Zhang .seealso: MatGetOwnershipRange(), MatMPISELLGetLocalMat() 1634d4002b98SHong Zhang 1635d4002b98SHong Zhang @*/ 1636d4002b98SHong Zhang PetscErrorCode MatMPISELLGetLocalMatCondensed(Mat A,MatReuse scall,IS *row,IS *col,Mat *A_loc) 1637d4002b98SHong Zhang { 1638d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 1639d4002b98SHong Zhang PetscErrorCode ierr; 1640d4002b98SHong Zhang PetscInt i,start,end,ncols,nzA,nzB,*cmap,imark,*idx; 1641d4002b98SHong Zhang IS isrowa,iscola; 1642d4002b98SHong Zhang Mat *aloc; 1643d4002b98SHong Zhang PetscBool match; 1644d4002b98SHong Zhang 1645d4002b98SHong Zhang PetscFunctionBegin; 1646d4002b98SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)A,MATMPISELL,&match);CHKERRQ(ierr); 1647d4002b98SHong Zhang if (!match) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP,"Requires MATMPISELL matrix as input"); 1648d4002b98SHong Zhang ierr = PetscLogEventBegin(MAT_Getlocalmatcondensed,A,0,0,0);CHKERRQ(ierr); 1649d4002b98SHong Zhang if (!row) { 1650d4002b98SHong Zhang start = A->rmap->rstart; end = A->rmap->rend; 1651d4002b98SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&isrowa);CHKERRQ(ierr); 1652d4002b98SHong Zhang } else { 1653d4002b98SHong Zhang isrowa = *row; 1654d4002b98SHong Zhang } 1655d4002b98SHong Zhang if (!col) { 1656d4002b98SHong Zhang start = A->cmap->rstart; 1657d4002b98SHong Zhang cmap = a->garray; 1658d4002b98SHong Zhang nzA = a->A->cmap->n; 1659d4002b98SHong Zhang nzB = a->B->cmap->n; 1660d4002b98SHong Zhang ierr = PetscMalloc1(nzA+nzB, &idx);CHKERRQ(ierr); 1661d4002b98SHong Zhang ncols = 0; 1662d4002b98SHong Zhang for (i=0; i<nzB; i++) { 1663d4002b98SHong Zhang if (cmap[i] < start) idx[ncols++] = cmap[i]; 1664d4002b98SHong Zhang else break; 1665d4002b98SHong Zhang } 1666d4002b98SHong Zhang imark = i; 1667d4002b98SHong Zhang for (i=0; i<nzA; i++) idx[ncols++] = start + i; 1668d4002b98SHong Zhang for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; 1669d4002b98SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&iscola);CHKERRQ(ierr); 1670d4002b98SHong Zhang } else { 1671d4002b98SHong Zhang iscola = *col; 1672d4002b98SHong Zhang } 1673d4002b98SHong Zhang if (scall != MAT_INITIAL_MATRIX) { 1674d4002b98SHong Zhang ierr = PetscMalloc1(1,&aloc);CHKERRQ(ierr); 1675d4002b98SHong Zhang aloc[0] = *A_loc; 1676d4002b98SHong Zhang } 1677d4002b98SHong Zhang ierr = MatCreateSubMatrices(A,1,&isrowa,&iscola,scall,&aloc);CHKERRQ(ierr); 1678d4002b98SHong Zhang *A_loc = aloc[0]; 1679d4002b98SHong Zhang ierr = PetscFree(aloc);CHKERRQ(ierr); 1680d4002b98SHong Zhang if (!row) { 1681d4002b98SHong Zhang ierr = ISDestroy(&isrowa);CHKERRQ(ierr); 1682d4002b98SHong Zhang } 1683d4002b98SHong Zhang if (!col) { 1684d4002b98SHong Zhang ierr = ISDestroy(&iscola);CHKERRQ(ierr); 1685d4002b98SHong Zhang } 1686d4002b98SHong Zhang ierr = PetscLogEventEnd(MAT_Getlocalmatcondensed,A,0,0,0);CHKERRQ(ierr); 1687d4002b98SHong Zhang PetscFunctionReturn(0); 1688d4002b98SHong Zhang } 1689d4002b98SHong Zhang 1690d4002b98SHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h> 1691d4002b98SHong Zhang 1692d4002b98SHong Zhang PetscErrorCode MatConvert_MPISELL_MPIAIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 1693d4002b98SHong Zhang { 1694d4002b98SHong Zhang PetscErrorCode ierr; 1695d4002b98SHong Zhang Mat_MPISELL *a=(Mat_MPISELL*)A->data; 1696d4002b98SHong Zhang Mat B; 1697d4002b98SHong Zhang Mat_MPIAIJ *b; 1698d4002b98SHong Zhang 1699d4002b98SHong Zhang PetscFunctionBegin; 1700d4002b98SHong Zhang if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Matrix must be assembled"); 1701d4002b98SHong Zhang 170294a8b381SRichard Tran Mills if (reuse == MAT_REUSE_MATRIX) { 170394a8b381SRichard Tran Mills B = *newmat; 170494a8b381SRichard Tran Mills } else { 1705d4002b98SHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 1706d4002b98SHong Zhang ierr = MatSetType(B,MATMPIAIJ);CHKERRQ(ierr); 1707d4002b98SHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 1708d4002b98SHong Zhang ierr = MatSetBlockSizes(B,A->rmap->bs,A->cmap->bs);CHKERRQ(ierr); 1709d4002b98SHong Zhang ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 1710d4002b98SHong Zhang ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 171194a8b381SRichard Tran Mills } 1712d4002b98SHong Zhang b = (Mat_MPIAIJ*) B->data; 171394a8b381SRichard Tran Mills 171494a8b381SRichard Tran Mills if (reuse == MAT_REUSE_MATRIX) { 171594a8b381SRichard Tran Mills ierr = MatConvert_SeqSELL_SeqAIJ(a->A, MATSEQAIJ, MAT_REUSE_MATRIX, &b->A);CHKERRQ(ierr); 171694a8b381SRichard Tran Mills ierr = MatConvert_SeqSELL_SeqAIJ(a->B, MATSEQAIJ, MAT_REUSE_MATRIX, &b->B);CHKERRQ(ierr); 171794a8b381SRichard Tran Mills } else { 1718d4002b98SHong Zhang ierr = MatDestroy(&b->A);CHKERRQ(ierr); 1719d4002b98SHong Zhang ierr = MatDestroy(&b->B);CHKERRQ(ierr); 1720d4002b98SHong Zhang ierr = MatDisAssemble_MPISELL(A);CHKERRQ(ierr); 1721d4002b98SHong Zhang ierr = MatConvert_SeqSELL_SeqAIJ(a->A, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->A);CHKERRQ(ierr); 1722d4002b98SHong Zhang ierr = MatConvert_SeqSELL_SeqAIJ(a->B, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->B);CHKERRQ(ierr); 1723d4002b98SHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1724d4002b98SHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1725d4002b98SHong Zhang ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1726d4002b98SHong Zhang ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 172794a8b381SRichard Tran Mills } 1728d4002b98SHong Zhang 1729d4002b98SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 1730d4002b98SHong Zhang ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr); 1731d4002b98SHong Zhang } else { 1732d4002b98SHong Zhang *newmat = B; 1733d4002b98SHong Zhang } 1734d4002b98SHong Zhang PetscFunctionReturn(0); 1735d4002b98SHong Zhang } 1736d4002b98SHong Zhang 1737d4002b98SHong Zhang PetscErrorCode MatConvert_MPIAIJ_MPISELL(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 1738d4002b98SHong Zhang { 1739d4002b98SHong Zhang PetscErrorCode ierr; 1740d4002b98SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 1741d4002b98SHong Zhang Mat B; 1742d4002b98SHong Zhang Mat_MPISELL *b; 1743d4002b98SHong Zhang 1744d4002b98SHong Zhang PetscFunctionBegin; 1745d4002b98SHong Zhang if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Matrix must be assembled"); 1746d4002b98SHong Zhang 174794a8b381SRichard Tran Mills if (reuse == MAT_REUSE_MATRIX) { 174894a8b381SRichard Tran Mills B = *newmat; 174994a8b381SRichard Tran Mills } else { 1750d4002b98SHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 1751d4002b98SHong Zhang ierr = MatSetType(B,MATMPISELL);CHKERRQ(ierr); 1752d4002b98SHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 1753d4002b98SHong Zhang ierr = MatSetBlockSizes(B,A->rmap->bs,A->cmap->bs);CHKERRQ(ierr); 1754d4002b98SHong Zhang ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 1755d4002b98SHong Zhang ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 175694a8b381SRichard Tran Mills } 1757d4002b98SHong Zhang b = (Mat_MPISELL*) B->data; 175894a8b381SRichard Tran Mills 175994a8b381SRichard Tran Mills if (reuse == MAT_REUSE_MATRIX) { 176094a8b381SRichard Tran Mills ierr = MatConvert_SeqAIJ_SeqSELL(a->A, MATSEQSELL, MAT_REUSE_MATRIX, &b->A);CHKERRQ(ierr); 176194a8b381SRichard Tran Mills ierr = MatConvert_SeqAIJ_SeqSELL(a->B, MATSEQSELL, MAT_REUSE_MATRIX, &b->B);CHKERRQ(ierr); 176294a8b381SRichard Tran Mills } else { 1763d4002b98SHong Zhang ierr = MatDestroy(&b->A);CHKERRQ(ierr); 1764d4002b98SHong Zhang ierr = MatDestroy(&b->B);CHKERRQ(ierr); 1765d4002b98SHong Zhang ierr = MatDisAssemble_MPIAIJ(A);CHKERRQ(ierr); 1766d4002b98SHong Zhang ierr = MatConvert_SeqAIJ_SeqSELL(a->A, MATSEQSELL, MAT_INITIAL_MATRIX, &b->A);CHKERRQ(ierr); 1767d4002b98SHong Zhang ierr = MatConvert_SeqAIJ_SeqSELL(a->B, MATSEQSELL, MAT_INITIAL_MATRIX, &b->B);CHKERRQ(ierr); 1768d4002b98SHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1769d4002b98SHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1770d4002b98SHong Zhang ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1771d4002b98SHong Zhang ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 177294a8b381SRichard Tran Mills } 1773d4002b98SHong Zhang 1774d4002b98SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 1775d4002b98SHong Zhang ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr); 1776d4002b98SHong Zhang } else { 1777d4002b98SHong Zhang *newmat = B; 1778d4002b98SHong Zhang } 1779d4002b98SHong Zhang PetscFunctionReturn(0); 1780d4002b98SHong Zhang } 1781d4002b98SHong Zhang 1782d4002b98SHong Zhang PetscErrorCode MatSOR_MPISELL(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 1783d4002b98SHong Zhang { 1784d4002b98SHong Zhang Mat_MPISELL *mat=(Mat_MPISELL*)matin->data; 1785d4002b98SHong Zhang PetscErrorCode ierr; 1786d4002b98SHong Zhang Vec bb1=0; 1787d4002b98SHong Zhang 1788d4002b98SHong Zhang PetscFunctionBegin; 1789d4002b98SHong Zhang if (flag == SOR_APPLY_UPPER) { 1790d4002b98SHong Zhang ierr = (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);CHKERRQ(ierr); 1791d4002b98SHong Zhang PetscFunctionReturn(0); 1792d4002b98SHong Zhang } 1793d4002b98SHong Zhang 1794d4002b98SHong Zhang if (its > 1 || ~flag & SOR_ZERO_INITIAL_GUESS || flag & SOR_EISENSTAT) { 1795d4002b98SHong Zhang ierr = VecDuplicate(bb,&bb1);CHKERRQ(ierr); 1796d4002b98SHong Zhang } 1797d4002b98SHong Zhang 1798d4002b98SHong Zhang if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) { 1799d4002b98SHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 1800d4002b98SHong Zhang ierr = (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);CHKERRQ(ierr); 1801d4002b98SHong Zhang its--; 1802d4002b98SHong Zhang } 1803d4002b98SHong Zhang 1804d4002b98SHong Zhang while (its--) { 1805d4002b98SHong Zhang ierr = VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1806d4002b98SHong Zhang ierr = VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1807d4002b98SHong Zhang 1808d4002b98SHong Zhang /* update rhs: bb1 = bb - B*x */ 1809d4002b98SHong Zhang ierr = VecScale(mat->lvec,-1.0);CHKERRQ(ierr); 1810d4002b98SHong Zhang ierr = (*mat->B->ops->multadd)(mat->B,mat->lvec,bb,bb1);CHKERRQ(ierr); 1811d4002b98SHong Zhang 1812d4002b98SHong Zhang /* local sweep */ 1813d4002b98SHong Zhang ierr = (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,1,xx);CHKERRQ(ierr); 1814d4002b98SHong Zhang } 1815d4002b98SHong Zhang } else if (flag & SOR_LOCAL_FORWARD_SWEEP) { 1816d4002b98SHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 1817d4002b98SHong Zhang ierr = (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);CHKERRQ(ierr); 1818d4002b98SHong Zhang its--; 1819d4002b98SHong Zhang } 1820d4002b98SHong Zhang while (its--) { 1821d4002b98SHong Zhang ierr = VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1822d4002b98SHong Zhang ierr = VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1823d4002b98SHong Zhang 1824d4002b98SHong Zhang /* update rhs: bb1 = bb - B*x */ 1825d4002b98SHong Zhang ierr = VecScale(mat->lvec,-1.0);CHKERRQ(ierr); 1826d4002b98SHong Zhang ierr = (*mat->B->ops->multadd)(mat->B,mat->lvec,bb,bb1);CHKERRQ(ierr); 1827d4002b98SHong Zhang 1828d4002b98SHong Zhang /* local sweep */ 1829d4002b98SHong Zhang ierr = (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_FORWARD_SWEEP,fshift,lits,1,xx);CHKERRQ(ierr); 1830d4002b98SHong Zhang } 1831d4002b98SHong Zhang } else if (flag & SOR_LOCAL_BACKWARD_SWEEP) { 1832d4002b98SHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 1833d4002b98SHong Zhang ierr = (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);CHKERRQ(ierr); 1834d4002b98SHong Zhang its--; 1835d4002b98SHong Zhang } 1836d4002b98SHong Zhang while (its--) { 1837d4002b98SHong Zhang ierr = VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1838d4002b98SHong Zhang ierr = VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1839d4002b98SHong Zhang 1840d4002b98SHong Zhang /* update rhs: bb1 = bb - B*x */ 1841d4002b98SHong Zhang ierr = VecScale(mat->lvec,-1.0);CHKERRQ(ierr); 1842d4002b98SHong Zhang ierr = (*mat->B->ops->multadd)(mat->B,mat->lvec,bb,bb1);CHKERRQ(ierr); 1843d4002b98SHong Zhang 1844d4002b98SHong Zhang /* local sweep */ 1845d4002b98SHong Zhang ierr = (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_BACKWARD_SWEEP,fshift,lits,1,xx);CHKERRQ(ierr); 1846d4002b98SHong Zhang } 1847d4002b98SHong Zhang } else SETERRQ(PetscObjectComm((PetscObject)matin),PETSC_ERR_SUP,"Parallel SOR not supported"); 1848d4002b98SHong Zhang 1849d4002b98SHong Zhang ierr = VecDestroy(&bb1);CHKERRQ(ierr); 1850d4002b98SHong Zhang 1851d4002b98SHong Zhang matin->factorerrortype = mat->A->factorerrortype; 1852d4002b98SHong Zhang PetscFunctionReturn(0); 1853d4002b98SHong Zhang } 1854d4002b98SHong Zhang 1855d4002b98SHong Zhang /*MC 1856d4002b98SHong Zhang MATMPISELL - MATMPISELL = "MPISELL" - A matrix type to be used for parallel sparse matrices. 1857d4002b98SHong Zhang 1858d4002b98SHong Zhang Options Database Keys: 1859d4002b98SHong Zhang . -mat_type MPISELL - sets the matrix type to "MPISELL" during a call to MatSetFromOptions() 1860d4002b98SHong Zhang 1861d4002b98SHong Zhang Level: beginner 1862d4002b98SHong Zhang 1863d4002b98SHong Zhang .seealso: MatCreateSELL() 1864d4002b98SHong Zhang M*/ 1865d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatCreate_MPISELL(Mat B) 1866d4002b98SHong Zhang { 1867d4002b98SHong Zhang Mat_MPISELL *b; 1868d4002b98SHong Zhang PetscErrorCode ierr; 1869d4002b98SHong Zhang PetscMPIInt size; 1870d4002b98SHong Zhang 1871d4002b98SHong Zhang PetscFunctionBegin; 1872d4002b98SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);CHKERRQ(ierr); 1873d4002b98SHong Zhang ierr = PetscNewLog(B,&b);CHKERRQ(ierr); 1874d4002b98SHong Zhang B->data = (void*)b; 1875d4002b98SHong Zhang ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1876d4002b98SHong Zhang B->assembled = PETSC_FALSE; 1877d4002b98SHong Zhang B->insertmode = NOT_SET_VALUES; 1878d4002b98SHong Zhang b->size = size; 1879d4002b98SHong Zhang ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)B),&b->rank);CHKERRQ(ierr); 1880d4002b98SHong Zhang /* build cache for off array entries formed */ 1881d4002b98SHong Zhang ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash);CHKERRQ(ierr); 1882d4002b98SHong Zhang 1883d4002b98SHong Zhang b->donotstash = PETSC_FALSE; 1884d4002b98SHong Zhang b->colmap = 0; 1885d4002b98SHong Zhang b->garray = 0; 1886d4002b98SHong Zhang b->roworiented = PETSC_TRUE; 1887d4002b98SHong Zhang 1888d4002b98SHong Zhang /* stuff used for matrix vector multiply */ 1889d4002b98SHong Zhang b->lvec = NULL; 1890d4002b98SHong Zhang b->Mvctx = NULL; 1891d4002b98SHong Zhang 1892d4002b98SHong Zhang /* stuff for MatGetRow() */ 1893d4002b98SHong Zhang b->rowindices = 0; 1894d4002b98SHong Zhang b->rowvalues = 0; 1895d4002b98SHong Zhang b->getrowactive = PETSC_FALSE; 1896d4002b98SHong Zhang 1897d4002b98SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_MPISELL);CHKERRQ(ierr); 1898d4002b98SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_MPISELL);CHKERRQ(ierr); 1899d4002b98SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatIsTranspose_C",MatIsTranspose_MPISELL);CHKERRQ(ierr); 1900d4002b98SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISELLSetPreallocation_C",MatMPISELLSetPreallocation_MPISELL);CHKERRQ(ierr); 1901d4002b98SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisell_mpiaij_C",MatConvert_MPISELL_MPIAIJ);CHKERRQ(ierr); 1902d4002b98SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatDiagonalScaleLocal_C",MatDiagonalScaleLocal_MPISELL);CHKERRQ(ierr); 1903d4002b98SHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPISELL);CHKERRQ(ierr); 1904d4002b98SHong Zhang PetscFunctionReturn(0); 1905d4002b98SHong Zhang } 1906