xref: /petsc/src/mat/impls/baij/seq/aijbaij.c (revision 7c307921e16164bcabd288ffaae1e907d3e04a89)
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"
8be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqBAIJ_SeqAIJ(Mat A,const 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;
13521d7252SBarry Smith   PetscInt       bs = A->bs,*ai = a->i,*aj = a->j,n = A->M/bs,i,j,k;
14690b6cddSBarry Smith   PetscInt       *rowlengths,*rows,*cols,maxlen = 0,ncols;
1587828ca2SBarry Smith   PetscScalar    *aa = a->a;
1661e22dc2SBarry Smith 
1761e22dc2SBarry Smith   PetscFunctionBegin;
18*7c307921SBarry 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   }
2518f87311SHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&B);CHKERRQ(ierr);
26f204ca49SKris Buschelman   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
27f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
28676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_COLUMN_ORIENTED);CHKERRQ(ierr);
29676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
30676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
3161e22dc2SBarry Smith   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
3261e22dc2SBarry Smith 
33*7c307921SBarry Smith   ierr = PetscMalloc(bs*sizeof(PetscInt),&rows);CHKERRQ(ierr);
34*7c307921SBarry Smith   ierr = PetscMalloc(bs*maxlen*sizeof(PetscInt),&cols);CHKERRQ(ierr);
35b9b97703SBarry Smith   for (i=0; i<n; i++) {
36b9b97703SBarry Smith     for (j=0; j<bs; j++) {
37b9b97703SBarry Smith       rows[j] = i*bs+j;
38b9b97703SBarry Smith     }
39b9b97703SBarry Smith     ncols = ai[i+1] - ai[i];
40b9b97703SBarry Smith     for (k=0; k<ncols; k++) {
41b9b97703SBarry Smith       for (j=0; j<bs; j++) {
42b9b97703SBarry Smith         cols[k*bs+j] = bs*(*aj) + j;
43b9b97703SBarry Smith       }
44b9b97703SBarry Smith       aj++;
45b9b97703SBarry Smith     }
46676c34cdSKris Buschelman     ierr  = MatSetValues(B,bs,rows,bs*ncols,cols,aa,INSERT_VALUES);CHKERRQ(ierr);
47b9b97703SBarry Smith     aa   += ncols*bs*bs;
48b9b97703SBarry Smith   }
49b9b97703SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
50b9b97703SBarry Smith   ierr = PetscFree(rows);CHKERRQ(ierr);
51676c34cdSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
52676c34cdSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
53676c34cdSKris Buschelman 
54521d7252SBarry Smith   B->bs = A->bs;
55521d7252SBarry Smith 
56ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
578ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
58c3d102feSKris Buschelman   } else {
5918f87311SHong Zhang     *newmat = B;
60c3d102feSKris Buschelman   }
6161e22dc2SBarry Smith   PetscFunctionReturn(0);
6261e22dc2SBarry Smith }
63273d9f13SBarry Smith EXTERN_C_END
6461e22dc2SBarry Smith 
6585fc7724SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
6685fc7724SBarry Smith 
6785fc7724SBarry Smith EXTERN_C_BEGIN
6885fc7724SBarry Smith #undef __FUNCT__
6985fc7724SBarry Smith #define __FUNCT__ "MatConvert_SeqAIJ_SeqBAIJ"
70be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat)
71521d7252SBarry Smith {
72676c34cdSKris Buschelman   Mat            B;
7385fc7724SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7485fc7724SBarry Smith   Mat_SeqBAIJ    *b;
75dfbe8321SBarry Smith   PetscErrorCode ierr;
76690b6cddSBarry Smith   PetscInt       *ai=a->i,m=A->M,n=A->N,i,*rowlengths;
7785fc7724SBarry Smith 
7885fc7724SBarry Smith   PetscFunctionBegin;
7985fc7724SBarry Smith   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
8085fc7724SBarry Smith 
81*7c307921SBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
8285fc7724SBarry Smith   for (i=0; i<m; i++) {
8385fc7724SBarry Smith     rowlengths[i] = ai[i+1] - ai[i];
8485fc7724SBarry Smith   }
8518f87311SHong Zhang   ierr = MatCreate(A->comm,m,n,m,n,&B);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 
90676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr);
91676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
92676c34cdSKris Buschelman   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
9385fc7724SBarry Smith 
94676c34cdSKris Buschelman   b  = (Mat_SeqBAIJ*)(B->data);
9585fc7724SBarry Smith 
96*7c307921SBarry Smith   ierr = PetscMemcpy(b->i,a->i,(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
97*7c307921SBarry Smith   ierr = PetscMemcpy(b->ilen,a->ilen,m*sizeof(PetscInt));CHKERRQ(ierr);
98*7c307921SBarry Smith   ierr = PetscMemcpy(b->j,a->j,a->nz*sizeof(PetscInt));CHKERRQ(ierr);
9985fc7724SBarry Smith   ierr = PetscMemcpy(b->a,a->a,a->nz*sizeof(MatScalar));CHKERRQ(ierr);
10085fc7724SBarry Smith 
101676c34cdSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
102676c34cdSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
103676c34cdSKris Buschelman 
104ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
1058ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
106c3d102feSKris Buschelman   } else {
10718f87311SHong Zhang     *newmat = B;
108c3d102feSKris Buschelman   }
10985fc7724SBarry Smith   PetscFunctionReturn(0);
11085fc7724SBarry Smith }
11185fc7724SBarry Smith EXTERN_C_END
112