1be1d678aSKris Buschelman #define PETSCMAT_DLL 261e22dc2SBarry Smith 361e22dc2SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 461e22dc2SBarry Smith 5273d9f13SBarry Smith EXTERN_C_BEGIN 64a2ae208SSatish Balay #undef __FUNCT__ 7a0aab4a4SKris Buschelman #define __FUNCT__ "MatConvert_SeqBAIJ_SeqAIJ" 8f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 9521d7252SBarry Smith { 10676c34cdSKris Buschelman Mat B; 1161e22dc2SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 12dfbe8321SBarry Smith PetscErrorCode ierr; 13*d0f46423SBarry Smith PetscInt bs = A->rmap->bs,*ai = a->i,*aj = a->j,n = A->rmap->N/bs,i,j,k; 14690b6cddSBarry Smith PetscInt *rowlengths,*rows,*cols,maxlen = 0,ncols; 15dd6ea824SBarry Smith MatScalar *aa = a->a; 1661e22dc2SBarry Smith 1761e22dc2SBarry Smith PetscFunctionBegin; 187c307921SBarry Smith ierr = PetscMalloc(n*bs*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 1961e22dc2SBarry Smith for (i=0; i<n; i++) { 20b9b97703SBarry Smith maxlen = PetscMax(maxlen,(ai[i+1] - ai[i])); 2161e22dc2SBarry Smith for (j=0; j<bs; j++) { 22b9b97703SBarry Smith rowlengths[i*bs+j] = bs*(ai[i+1] - ai[i]); 2361e22dc2SBarry Smith } 2461e22dc2SBarry Smith } 257adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 26*d0f46423SBarry Smith ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 27f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 28f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr); 294e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 3061e22dc2SBarry Smith ierr = PetscFree(rowlengths);CHKERRQ(ierr); 3161e22dc2SBarry Smith 327c307921SBarry Smith ierr = PetscMalloc(bs*sizeof(PetscInt),&rows);CHKERRQ(ierr); 337c307921SBarry Smith ierr = PetscMalloc(bs*maxlen*sizeof(PetscInt),&cols);CHKERRQ(ierr); 34b9b97703SBarry Smith for (i=0; i<n; i++) { 35b9b97703SBarry Smith for (j=0; j<bs; j++) { 36b9b97703SBarry Smith rows[j] = i*bs+j; 37b9b97703SBarry Smith } 38b9b97703SBarry Smith ncols = ai[i+1] - ai[i]; 39b9b97703SBarry Smith for (k=0; k<ncols; k++) { 40b9b97703SBarry Smith for (j=0; j<bs; j++) { 41b9b97703SBarry Smith cols[k*bs+j] = bs*(*aj) + j; 42b9b97703SBarry Smith } 43b9b97703SBarry Smith aj++; 44b9b97703SBarry Smith } 45676c34cdSKris Buschelman ierr = MatSetValues(B,bs,rows,bs*ncols,cols,aa,INSERT_VALUES);CHKERRQ(ierr); 46b9b97703SBarry Smith aa += ncols*bs*bs; 47b9b97703SBarry Smith } 48b9b97703SBarry Smith ierr = PetscFree(cols);CHKERRQ(ierr); 49b9b97703SBarry Smith ierr = PetscFree(rows);CHKERRQ(ierr); 50676c34cdSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 51676c34cdSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 52676c34cdSKris Buschelman 53*d0f46423SBarry Smith B->rmap->bs = A->rmap->bs; 54521d7252SBarry Smith 55ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 568ab5b326SKris Buschelman ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 57c3d102feSKris Buschelman } else { 5818f87311SHong Zhang *newmat = B; 59c3d102feSKris Buschelman } 6061e22dc2SBarry Smith PetscFunctionReturn(0); 6161e22dc2SBarry Smith } 62273d9f13SBarry Smith EXTERN_C_END 6361e22dc2SBarry Smith 6485fc7724SBarry Smith #include "src/mat/impls/aij/seq/aij.h" 6585fc7724SBarry Smith 6685fc7724SBarry Smith EXTERN_C_BEGIN 6785fc7724SBarry Smith #undef __FUNCT__ 6885fc7724SBarry Smith #define __FUNCT__ "MatConvert_SeqAIJ_SeqBAIJ" 698cf70c4bSSatish Balay PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) 70521d7252SBarry Smith { 71676c34cdSKris Buschelman Mat B; 7285fc7724SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7385fc7724SBarry Smith Mat_SeqBAIJ *b; 74dfbe8321SBarry Smith PetscErrorCode ierr; 75*d0f46423SBarry Smith PetscInt *ai=a->i,m=A->rmap->N,n=A->cmap->N,i,*rowlengths; 7685fc7724SBarry Smith 7785fc7724SBarry Smith PetscFunctionBegin; 7885fc7724SBarry Smith if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square"); 7985fc7724SBarry Smith 807c307921SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 8185fc7724SBarry Smith for (i=0; i<m; i++) { 8285fc7724SBarry Smith rowlengths[i] = ai[i+1] - ai[i]; 8385fc7724SBarry Smith } 847adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 85f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 86f204ca49SKris Buschelman ierr = MatSetType(B,newtype);CHKERRQ(ierr); 87ab93d7beSBarry Smith ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,1,0,rowlengths);CHKERRQ(ierr); 8885fc7724SBarry Smith ierr = PetscFree(rowlengths);CHKERRQ(ierr); 8985fc7724SBarry Smith 904e0d8c25SBarry Smith ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 9185fc7724SBarry Smith 92676c34cdSKris Buschelman b = (Mat_SeqBAIJ*)(B->data); 9385fc7724SBarry Smith 947c307921SBarry Smith ierr = PetscMemcpy(b->i,a->i,(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 957c307921SBarry Smith ierr = PetscMemcpy(b->ilen,a->ilen,m*sizeof(PetscInt));CHKERRQ(ierr); 967c307921SBarry Smith ierr = PetscMemcpy(b->j,a->j,a->nz*sizeof(PetscInt));CHKERRQ(ierr); 9785fc7724SBarry Smith ierr = PetscMemcpy(b->a,a->a,a->nz*sizeof(MatScalar));CHKERRQ(ierr); 9885fc7724SBarry Smith 99676c34cdSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 100676c34cdSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 101676c34cdSKris Buschelman 102ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 1038ab5b326SKris Buschelman ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 104c3d102feSKris Buschelman } else { 10518f87311SHong Zhang *newmat = B; 106c3d102feSKris Buschelman } 10785fc7724SBarry Smith PetscFunctionReturn(0); 10885fc7724SBarry Smith } 10985fc7724SBarry Smith EXTERN_C_END 110