xref: /petsc/src/mat/impls/sbaij/mpi/mpiaijsbaij.c (revision ad5247fdeeae2f8b9416c4cfe38baedadecec029)
1 
2 #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> /*I "petscmat.h" I*/
3 #include <../src/mat/impls/aij/mpi/mpiaij.h>
4 #include <private/matimpl.h>
5 #include <petscmat.h>
6 
7 EXTERN_C_BEGIN
8 #undef __FUNCT__
9 #define __FUNCT__ "MatConvert_MPIAIJ_MPISBAIJ"
10 PetscErrorCode  MatConvert_MPIAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
11 {
12   PetscErrorCode     ierr;
13   Mat                M;
14   Mat_MPIAIJ         *mpimat = (Mat_MPIAIJ*)A->data;
15   Mat_SeqAIJ         *Aa = (Mat_SeqAIJ*)mpimat->A->data,*Ba = (Mat_SeqAIJ*)mpimat->B->data;
16   PetscInt           *d_nnz,*o_nnz;
17   PetscInt           i,j,nz;
18   PetscInt           m,n,lm,ln;
19   PetscInt           rstart,rend;
20   const PetscScalar  *vwork;
21   const PetscInt     *cwork;
22 
23   PetscFunctionBegin;
24   if (!A->symmetric) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
25   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
26   ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr);
27   ierr = PetscMalloc2(lm,PetscInt,&d_nnz,lm,PetscInt,&o_nnz);CHKERRQ(ierr);
28 
29   ierr = MatMarkDiagonal_SeqAIJ(mpimat->A);CHKERRQ(ierr);
30   for(i=0;i<lm;i++){
31     d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
32     o_nnz[i] = Ba->i[i+1] - Ba->i[i];
33   }
34 
35   ierr = MatCreate(((PetscObject)A)->comm,&M);CHKERRQ(ierr);
36   ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
37   ierr = MatSetType(M,newtype);CHKERRQ(ierr);
38   ierr = MatMPISBAIJSetPreallocation(M,1,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
39 
40   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
41 
42   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
43   for(i=rstart;i<rend;i++){
44     ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
45     j = 0;
46     while (cwork[j] < i){ j++; nz--;}
47     ierr = MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);CHKERRQ(ierr);
48     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
49   }
50   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
51   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
52 
53   if (A->hermitian){
54     ierr = MatSetOption(M,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr);
55   }
56 
57   if (reuse == MAT_REUSE_MATRIX) {
58     ierr = MatHeaderReplace(A,M);CHKERRQ(ierr);
59   } else {
60     *newmat = M;
61   }
62   PetscFunctionReturn(0);
63 }
64 /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */
65 #undef __FUNCT__
66 #define __FUNCT__ "MatConvert_MPIBAIJ_MPISBAIJ"
67 PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
68 {
69   PetscErrorCode     ierr;
70   Mat                M;
71   Mat_MPIBAIJ        *mpimat = (Mat_MPIBAIJ*)A->data;
72   Mat_SeqBAIJ        *Aa = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data;
73   PetscInt           *d_nnz,*o_nnz;
74   PetscInt           i,j,nz;
75   PetscInt           m,n,lm,ln;
76   PetscInt           rstart,rend;
77   const PetscScalar  *vwork;
78   const PetscInt     *cwork;
79   PetscInt           bs = A->rmap->bs;
80 
81   PetscFunctionBegin;
82   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
83   ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr);
84   printf(" --- in MatConvert_MPIBAIJ_MPISBAIJ, bs = %d m %d lm %d\n", bs,m,lm);
85   ierr = PetscMalloc2(lm/bs,PetscInt,&d_nnz,lm/bs,PetscInt,&o_nnz);CHKERRQ(ierr);
86 
87   ierr = MatMarkDiagonal_SeqBAIJ(mpimat->A);CHKERRQ(ierr);
88   for(i=0;i<lm/bs;i++){
89     d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
90     o_nnz[i] = Ba->i[i+1] - Ba->i[i];
91   }
92 
93   ierr = MatCreate(((PetscObject)A)->comm,&M);CHKERRQ(ierr);
94   ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
95   ierr = MatSetType(M,newtype);CHKERRQ(ierr);
96   ierr = MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
97 
98   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
99 
100   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
101   printf("--- rstart  %d end %d\n",rstart,rend);
102   for(i=rstart;i<rend;i++){
103     ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
104     j = 0;
105     while (cwork[j] < i){ j++; nz--;}
106     ierr = MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);CHKERRQ(ierr); // missing lower triangular entries in the diagonal blocks!!!
107     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
108   }
109   printf(" --- ok 2\n");
110   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
111   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
112 
113   if (A->hermitian){
114     ierr = MatSetOption(M,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr);
115   }
116 
117   printf(" --- ok 3\n");
118   if (reuse == MAT_REUSE_MATRIX) {
119     ierr = MatHeaderReplace(A,M);CHKERRQ(ierr);
120   } else {
121     *newmat = M;
122   }
123   PetscFunctionReturn(0);
124 }
125 EXTERN_C_END
126