xref: /petsc/src/mat/impls/baij/seq/aijbaij.c (revision 2fc52814d27bf1f4e71021c1c3ebb532b583ed60)
1 /*$Id: aijbaij.c,v 1.9 2001/08/07 03:02:55 balay Exp $*/
2 
3 #include "src/mat/impls/baij/seq/baij.h"
4 
5 EXTERN_C_BEGIN
6 #undef __FUNCT__
7 #define __FUNCT__ "MatConvert_SeqBAI_SeqAIJ"
8 int MatConvert_SeqBAIJ_SeqAIJ(Mat A,const MatType newtype,Mat *newmat) {
9   Mat          B;
10   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
11   int          ierr,bs = a->bs,*ai = a->i,*aj = a->j,n = A->M/bs,i,j,k;
12   int          *rowlengths,*rows,*cols,maxlen = 0,ncols;
13   PetscScalar  *aa = a->a;
14 
15   PetscFunctionBegin;
16   ierr = PetscMalloc(n*bs*sizeof(int),&rowlengths);CHKERRQ(ierr);
17   for (i=0; i<n; i++) {
18     maxlen = PetscMax(maxlen,(ai[i+1] - ai[i]));
19     for (j=0; j<bs; j++) {
20       rowlengths[i*bs+j] = bs*(ai[i+1] - ai[i]);
21     }
22   }
23   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,A->m,A->n,0,rowlengths,&B);CHKERRQ(ierr);
24   ierr = MatSetOption(B,MAT_COLUMN_ORIENTED);CHKERRQ(ierr);
25   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
26   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
27   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
28 
29   ierr = PetscMalloc(bs*sizeof(int),&rows);CHKERRQ(ierr);
30   ierr = PetscMalloc(bs*maxlen*sizeof(int),&cols);CHKERRQ(ierr);
31   for (i=0; i<n; i++) {
32     for (j=0; j<bs; j++) {
33       rows[j] = i*bs+j;
34     }
35     ncols = ai[i+1] - ai[i];
36     for (k=0; k<ncols; k++) {
37       for (j=0; j<bs; j++) {
38         cols[k*bs+j] = bs*(*aj) + j;
39       }
40       aj++;
41     }
42     ierr  = MatSetValues(B,bs,rows,bs*ncols,cols,aa,INSERT_VALUES);CHKERRQ(ierr);
43     aa   += ncols*bs*bs;
44   }
45   ierr = PetscFree(cols);CHKERRQ(ierr);
46   ierr = PetscFree(rows);CHKERRQ(ierr);
47   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
49 
50   /* Fake support for "inplace" convert. */
51   if (*newmat == A) {
52     ierr = MatDestroy(A);CHKERRQ(ierr);
53   }
54   *newmat = B;
55   PetscFunctionReturn(0);
56 }
57 EXTERN_C_END
58 
59 #include "src/mat/impls/aij/seq/aij.h"
60 
61 EXTERN_C_BEGIN
62 #undef __FUNCT__
63 #define __FUNCT__ "MatConvert_SeqAIJ_SeqBAIJ"
64 int MatConvert_SeqAIJ_SeqBAIJ(Mat A,const MatType newtype,Mat *newmat) {
65   Mat         B;
66   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)A->data;
67   Mat_SeqBAIJ *b;
68   int         ierr,*ai=a->i,m=A->M,n=A->N,i,*rowlengths;
69 
70   PetscFunctionBegin;
71   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
72 
73   ierr = PetscMalloc(m*sizeof(int),&rowlengths);CHKERRQ(ierr);
74   for (i=0; i<m; i++) {
75     rowlengths[i] = ai[i+1] - ai[i];
76   }
77   ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,1,m,n,0,rowlengths,&B);CHKERRQ(ierr);
78   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
79 
80   ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr);
81   ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr);
82   ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
83 
84   b  = (Mat_SeqBAIJ*)(B->data);
85 
86   ierr = PetscMemcpy(b->i,a->i,(m+1)*sizeof(int));CHKERRQ(ierr);
87   ierr = PetscMemcpy(b->ilen,a->ilen,m*sizeof(int));CHKERRQ(ierr);
88   ierr = PetscMemcpy(b->j,a->j,a->nz*sizeof(int));CHKERRQ(ierr);
89   ierr = PetscMemcpy(b->a,a->a,a->nz*sizeof(MatScalar));CHKERRQ(ierr);
90 
91   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
93 
94   /* Fake support for "inplace" convert. */
95   if (*newmat == A) {
96     ierr = MatDestroy(A);CHKERRQ(ierr);
97   }
98   *newmat = B;
99   PetscFunctionReturn(0);
100 }
101 EXTERN_C_END
102