161e22dc2SBarry Smith 2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> 361e22dc2SBarry Smith 44a2ae208SSatish Balay #undef __FUNCT__ 5a0aab4a4SKris Buschelman #define __FUNCT__ "MatConvert_SeqBAIJ_SeqAIJ" 68cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 7521d7252SBarry Smith { 8676c34cdSKris Buschelman Mat B; 961e22dc2SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 10dfbe8321SBarry Smith PetscErrorCode ierr; 11d0f46423SBarry Smith PetscInt bs = A->rmap->bs,*ai = a->i,*aj = a->j,n = A->rmap->N/bs,i,j,k; 12690b6cddSBarry Smith PetscInt *rowlengths,*rows,*cols,maxlen = 0,ncols; 13dd6ea824SBarry Smith MatScalar *aa = a->a; 1461e22dc2SBarry Smith 1561e22dc2SBarry Smith PetscFunctionBegin; 16*785e854fSJed Brown ierr = PetscMalloc1(n*bs,&rowlengths);CHKERRQ(ierr); 1761e22dc2SBarry Smith for (i=0; i<n; i++) { 18b9b97703SBarry Smith maxlen = PetscMax(maxlen,(ai[i+1] - ai[i])); 1961e22dc2SBarry Smith for (j=0; j<bs; j++) { 20b9b97703SBarry Smith rowlengths[i*bs+j] = bs*(ai[i+1] - ai[i]); 2161e22dc2SBarry Smith } 2261e22dc2SBarry Smith } 23ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 24d0f46423SBarry Smith ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 256d0a4a0eSHong Zhang ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr); 26f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr); 274e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 2861e22dc2SBarry Smith ierr = PetscFree(rowlengths);CHKERRQ(ierr); 2961e22dc2SBarry Smith 30*785e854fSJed Brown ierr = PetscMalloc1(bs,&rows);CHKERRQ(ierr); 31*785e854fSJed Brown ierr = PetscMalloc1(bs*maxlen,&cols);CHKERRQ(ierr); 32b9b97703SBarry Smith for (i=0; i<n; i++) { 33b9b97703SBarry Smith for (j=0; j<bs; j++) { 34b9b97703SBarry Smith rows[j] = i*bs+j; 35b9b97703SBarry Smith } 36b9b97703SBarry Smith ncols = ai[i+1] - ai[i]; 37b9b97703SBarry Smith for (k=0; k<ncols; k++) { 38b9b97703SBarry Smith for (j=0; j<bs; j++) { 39b9b97703SBarry Smith cols[k*bs+j] = bs*(*aj) + j; 40b9b97703SBarry Smith } 41b9b97703SBarry Smith aj++; 42b9b97703SBarry Smith } 43676c34cdSKris Buschelman ierr = MatSetValues(B,bs,rows,bs*ncols,cols,aa,INSERT_VALUES);CHKERRQ(ierr); 44b9b97703SBarry Smith aa += ncols*bs*bs; 45b9b97703SBarry Smith } 46b9b97703SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 47b9b97703SBarry Smith ierr = PetscFree(rows);CHKERRQ(ierr); 48676c34cdSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 49676c34cdSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 50676c34cdSKris Buschelman 51d0f46423SBarry Smith B->rmap->bs = A->rmap->bs; 52521d7252SBarry Smith 53ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 548ab5b326SKris Buschelman ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 55c3d102feSKris Buschelman } else { 5618f87311SHong Zhang *newmat = B; 57c3d102feSKris Buschelman } 5861e22dc2SBarry Smith PetscFunctionReturn(0); 5961e22dc2SBarry Smith } 6061e22dc2SBarry Smith 61c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 6285fc7724SBarry Smith 6385fc7724SBarry Smith #undef __FUNCT__ 6485fc7724SBarry Smith #define __FUNCT__ "MatConvert_SeqAIJ_SeqBAIJ" 658cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 66521d7252SBarry Smith { 67676c34cdSKris Buschelman Mat B; 6885fc7724SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 6985fc7724SBarry Smith Mat_SeqBAIJ *b; 70dfbe8321SBarry Smith PetscErrorCode ierr; 71d0f46423SBarry Smith PetscInt *ai=a->i,m=A->rmap->N,n=A->cmap->N,i,*rowlengths; 7285fc7724SBarry Smith 7385fc7724SBarry Smith PetscFunctionBegin; 74e32f2f54SBarry Smith if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square"); 7585fc7724SBarry Smith 76*785e854fSJed Brown ierr = PetscMalloc1(m,&rowlengths);CHKERRQ(ierr); 7785fc7724SBarry Smith for (i=0; i<m; i++) { 7885fc7724SBarry Smith rowlengths[i] = ai[i+1] - ai[i]; 7985fc7724SBarry Smith } 80ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 81f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 826d0a4a0eSHong Zhang ierr = MatSetType(B,MATSEQBAIJ);CHKERRQ(ierr); 83ab93d7beSBarry Smith ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,1,0,rowlengths);CHKERRQ(ierr); 8485fc7724SBarry Smith ierr = PetscFree(rowlengths);CHKERRQ(ierr); 8585fc7724SBarry Smith 864e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 8785fc7724SBarry Smith 88676c34cdSKris Buschelman b = (Mat_SeqBAIJ*)(B->data); 8985fc7724SBarry Smith 907c307921SBarry Smith ierr = PetscMemcpy(b->i,a->i,(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 917c307921SBarry Smith ierr = PetscMemcpy(b->ilen,a->ilen,m*sizeof(PetscInt));CHKERRQ(ierr); 927c307921SBarry Smith ierr = PetscMemcpy(b->j,a->j,a->nz*sizeof(PetscInt));CHKERRQ(ierr); 9385fc7724SBarry Smith ierr = PetscMemcpy(b->a,a->a,a->nz*sizeof(MatScalar));CHKERRQ(ierr); 9485fc7724SBarry Smith 95676c34cdSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 96676c34cdSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 97676c34cdSKris Buschelman 98ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 998ab5b326SKris Buschelman ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 100c3d102feSKris Buschelman } else { 10118f87311SHong Zhang *newmat = B; 102c3d102feSKris Buschelman } 10385fc7724SBarry Smith PetscFunctionReturn(0); 10485fc7724SBarry Smith } 105