xref: /petsc/src/mat/impls/baij/seq/aijbaij.c (revision 85fc7724dce558a038329fe25788132d65e2c832)
173f4d377SMatthew Knepley /*$Id: aijbaij.c,v 1.9 2001/08/07 03:02:55 balay Exp $*/
261e22dc2SBarry Smith 
361e22dc2SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
461e22dc2SBarry Smith 
5273d9f13SBarry Smith EXTERN_C_BEGIN
64a2ae208SSatish Balay #undef __FUNCT__
74a2ae208SSatish Balay #define __FUNCT__ "MatConvert_SeqBAI_SeqAIJ"
8b9b97703SBarry Smith int MatConvert_SeqBAIJ_SeqAIJ(Mat A,MatType newtype,Mat *B)
961e22dc2SBarry Smith {
1061e22dc2SBarry Smith   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
11b9b97703SBarry Smith   int          ierr,bs = a->bs,*ai = a->i,*aj = a->j,n = A->M/bs,i,j,k;
12b9b97703SBarry Smith   int          *rowlengths,*rows,*cols,maxlen = 0,ncols;
1387828ca2SBarry Smith   PetscScalar  *aa = a->a;
1461e22dc2SBarry Smith 
1561e22dc2SBarry Smith   PetscFunctionBegin;
1682502324SSatish Balay   ierr = PetscMalloc(n*bs*sizeof(int),&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   }
23b9b97703SBarry Smith   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,A->m,A->n,0,rowlengths,B);CHKERRQ(ierr);
24b9b97703SBarry Smith   ierr = MatSetOption(*B,MAT_COLUMN_ORIENTED);CHKERRQ(ierr);
25b9b97703SBarry Smith   ierr = MatSetOption(*B,MAT_ROWS_SORTED);CHKERRQ(ierr);
26b9b97703SBarry Smith   ierr = MatSetOption(*B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
2761e22dc2SBarry Smith   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
2861e22dc2SBarry Smith 
2982502324SSatish Balay   ierr = PetscMalloc(bs*sizeof(int),&rows);CHKERRQ(ierr);
3082502324SSatish Balay   ierr = PetscMalloc(bs*maxlen*sizeof(int),&cols);CHKERRQ(ierr);
31b9b97703SBarry Smith   for (i=0; i<n; i++) {
32b9b97703SBarry Smith     for (j=0; j<bs; j++) {
33b9b97703SBarry Smith       rows[j] = i*bs+j;
34b9b97703SBarry Smith     }
35b9b97703SBarry Smith     ncols = ai[i+1] - ai[i];
36b9b97703SBarry Smith     for (k=0; k<ncols; k++) {
37b9b97703SBarry Smith       for (j=0; j<bs; j++) {
38b9b97703SBarry Smith         cols[k*bs+j] = bs*(*aj) + j;
39b9b97703SBarry Smith       }
40b9b97703SBarry Smith       aj++;
41b9b97703SBarry Smith     }
421fced5e8SSatish Balay     ierr  = MatSetValues(*B,bs,rows,bs*ncols,cols,aa,INSERT_VALUES);CHKERRQ(ierr);
43b9b97703SBarry Smith     aa   += ncols*bs*bs;
44b9b97703SBarry Smith   }
45b9b97703SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
46b9b97703SBarry Smith   ierr = PetscFree(rows);CHKERRQ(ierr);
47b9b97703SBarry Smith   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48b9b97703SBarry Smith   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4961e22dc2SBarry Smith   PetscFunctionReturn(0);
5061e22dc2SBarry Smith }
51273d9f13SBarry Smith EXTERN_C_END
5261e22dc2SBarry Smith 
53*85fc7724SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
54*85fc7724SBarry Smith 
55*85fc7724SBarry Smith EXTERN_C_BEGIN
56*85fc7724SBarry Smith #undef __FUNCT__
57*85fc7724SBarry Smith #define __FUNCT__ "MatConvert_SeqAIJ_SeqBAIJ"
58*85fc7724SBarry Smith int MatConvert_SeqAIJ_SeqBAIJ(Mat A,MatType newtype,Mat *B)
59*85fc7724SBarry Smith {
60*85fc7724SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)A->data;
61*85fc7724SBarry Smith   Mat_SeqBAIJ *b;
62*85fc7724SBarry Smith   int         ierr,*ai=a->i,m=A->M,n=A->N,i,*rowlengths;
63*85fc7724SBarry Smith 
64*85fc7724SBarry Smith   PetscFunctionBegin;
65*85fc7724SBarry Smith   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
66*85fc7724SBarry Smith 
67*85fc7724SBarry Smith   ierr = PetscMalloc(m*sizeof(int),&rowlengths);CHKERRQ(ierr);
68*85fc7724SBarry Smith   for (i=0; i<m; i++) {
69*85fc7724SBarry Smith     rowlengths[i] = ai[i+1] - ai[i];
70*85fc7724SBarry Smith   }
71*85fc7724SBarry Smith   ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,1,m,n,0,rowlengths,B);CHKERRQ(ierr);
72*85fc7724SBarry Smith   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
73*85fc7724SBarry Smith 
74*85fc7724SBarry Smith   ierr = MatSetOption(*B,MAT_ROW_ORIENTED);CHKERRQ(ierr);
75*85fc7724SBarry Smith   ierr = MatSetOption(*B,MAT_ROWS_SORTED);CHKERRQ(ierr);
76*85fc7724SBarry Smith   ierr = MatSetOption(*B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
77*85fc7724SBarry Smith 
78*85fc7724SBarry Smith   b  = (Mat_SeqBAIJ*)(*B)->data;
79*85fc7724SBarry Smith 
80*85fc7724SBarry Smith   ierr = PetscMemcpy(b->i,a->i,(m+1)*sizeof(int));CHKERRQ(ierr);
81*85fc7724SBarry Smith   ierr = PetscMemcpy(b->ilen,a->ilen,m*sizeof(int));CHKERRQ(ierr);
82*85fc7724SBarry Smith   ierr = PetscMemcpy(b->j,a->j,a->nz*sizeof(int));CHKERRQ(ierr);
83*85fc7724SBarry Smith   ierr = PetscMemcpy(b->a,a->a,a->nz*sizeof(MatScalar));CHKERRQ(ierr);
84*85fc7724SBarry Smith 
85*85fc7724SBarry Smith   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
86*85fc7724SBarry Smith   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
87*85fc7724SBarry Smith   PetscFunctionReturn(0);
88*85fc7724SBarry Smith }
89*85fc7724SBarry Smith EXTERN_C_END
90