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