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