xref: /petsc/src/mat/impls/baij/seq/aijbaij.c (revision 18f873111d4b6eb34a7a3679fe2e9fdf652c10bd)
161e22dc2SBarry Smith 
261e22dc2SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
361e22dc2SBarry Smith 
4273d9f13SBarry Smith EXTERN_C_BEGIN
54a2ae208SSatish Balay #undef __FUNCT__
64a2ae208SSatish Balay #define __FUNCT__ "MatConvert_SeqBAI_SeqAIJ"
7521d7252SBarry Smith PetscErrorCode MatConvert_SeqBAIJ_SeqAIJ(Mat A,const MatType newtype,Mat *newmat)
8521d7252SBarry Smith {
9676c34cdSKris Buschelman   Mat            B;
1061e22dc2SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
11dfbe8321SBarry Smith   PetscErrorCode ierr;
12521d7252SBarry Smith   PetscInt       bs = A->bs,*ai = a->i,*aj = a->j,n = A->M/bs,i,j,k;
13690b6cddSBarry Smith   PetscInt       *rowlengths,*rows,*cols,maxlen = 0,ncols;
1487828ca2SBarry Smith   PetscScalar    *aa = a->a;
1561e22dc2SBarry Smith 
1661e22dc2SBarry Smith   PetscFunctionBegin;
1782502324SSatish Balay   ierr = PetscMalloc(n*bs*sizeof(int),&rowlengths);CHKERRQ(ierr);
1861e22dc2SBarry Smith   for (i=0; i<n; i++) {
19b9b97703SBarry Smith     maxlen = PetscMax(maxlen,(ai[i+1] - ai[i]));
2061e22dc2SBarry Smith     for (j=0; j<bs; j++) {
21b9b97703SBarry Smith       rowlengths[i*bs+j] = bs*(ai[i+1] - ai[i]);
2261e22dc2SBarry Smith     }
2361e22dc2SBarry Smith   }
24*18f87311SHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&B);CHKERRQ(ierr);
25f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
26f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
27676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_COLUMN_ORIENTED);CHKERRQ(ierr);
28676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
29676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
3061e22dc2SBarry Smith   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
3161e22dc2SBarry Smith 
3282502324SSatish Balay   ierr = PetscMalloc(bs*sizeof(int),&rows);CHKERRQ(ierr);
3382502324SSatish Balay   ierr = PetscMalloc(bs*maxlen*sizeof(int),&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 
53521d7252SBarry Smith   B->bs = A->bs;
54521d7252SBarry Smith 
55676c34cdSKris Buschelman   /* Fake support for "inplace" convert. */
56*18f87311SHong Zhang   if (*newmat == A) {
57676c34cdSKris Buschelman     ierr = MatDestroy(A);CHKERRQ(ierr);
58676c34cdSKris Buschelman   }
59*18f87311SHong Zhang   *newmat = B;
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"
69521d7252SBarry Smith PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat A,const MatType newtype,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;
75690b6cddSBarry Smith   PetscInt       *ai=a->i,m=A->M,n=A->N,i,*rowlengths;
7685fc7724SBarry Smith 
7785fc7724SBarry Smith   PetscFunctionBegin;
7885fc7724SBarry Smith   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
7985fc7724SBarry Smith 
8085fc7724SBarry Smith   ierr = PetscMalloc(m*sizeof(int),&rowlengths);CHKERRQ(ierr);
8185fc7724SBarry Smith   for (i=0; i<m; i++) {
8285fc7724SBarry Smith     rowlengths[i] = ai[i+1] - ai[i];
8385fc7724SBarry Smith   }
84*18f87311SHong Zhang   ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr);
85f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
86f204ca49SKris Buschelman   ierr = MatSeqBAIJSetPreallocation(B,1,0,rowlengths);CHKERRQ(ierr);
8785fc7724SBarry Smith   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
8885fc7724SBarry Smith 
89676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr);
90676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
91676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
9285fc7724SBarry Smith 
93676c34cdSKris Buschelman   b  = (Mat_SeqBAIJ*)(B->data);
9485fc7724SBarry Smith 
9585fc7724SBarry Smith   ierr = PetscMemcpy(b->i,a->i,(m+1)*sizeof(int));CHKERRQ(ierr);
9685fc7724SBarry Smith   ierr = PetscMemcpy(b->ilen,a->ilen,m*sizeof(int));CHKERRQ(ierr);
9785fc7724SBarry Smith   ierr = PetscMemcpy(b->j,a->j,a->nz*sizeof(int));CHKERRQ(ierr);
9885fc7724SBarry Smith   ierr = PetscMemcpy(b->a,a->a,a->nz*sizeof(MatScalar));CHKERRQ(ierr);
9985fc7724SBarry Smith 
100676c34cdSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101676c34cdSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
102676c34cdSKris Buschelman 
103676c34cdSKris Buschelman   /* Fake support for "inplace" convert. */
104*18f87311SHong Zhang   if (*newmat == A) {
105676c34cdSKris Buschelman     ierr = MatDestroy(A);CHKERRQ(ierr);
106676c34cdSKris Buschelman   }
107*18f87311SHong Zhang   *newmat = B;
10885fc7724SBarry Smith   PetscFunctionReturn(0);
10985fc7724SBarry Smith }
11085fc7724SBarry Smith EXTERN_C_END
111