xref: /petsc/src/mat/impls/sbaij/mpi/mpiaijsbaij.c (revision e27a552be151e08936ff7d65f1f2e8dae4b63b83)
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 (reuse == MAT_REUSE_MATRIX) {
54     ierr = MatHeaderReplace(A,M);CHKERRQ(ierr);
55   } else {
56     *newmat = M;
57   }
58   PetscFunctionReturn(0);
59 }
60 /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */
61 #undef __FUNCT__
62 #define __FUNCT__ "MatConvert_MPIBAIJ_MPISBAIJ"
63 PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
64 {
65   PetscErrorCode     ierr;
66   Mat                M;
67   Mat_MPIBAIJ        *mpimat = (Mat_MPIBAIJ*)A->data;
68   Mat_SeqBAIJ        *Aa = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data;
69   PetscInt           *d_nnz,*o_nnz;
70   PetscInt           i,j,nz;
71   PetscInt           m,n,lm,ln;
72   PetscInt           rstart,rend;
73   const PetscScalar  *vwork;
74   const PetscInt     *cwork;
75   PetscInt           bs = A->rmap->bs;
76 
77   PetscFunctionBegin;
78   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
79   ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr);
80   ierr = PetscMalloc2(lm/bs,PetscInt,&d_nnz,lm/bs,PetscInt,&o_nnz);CHKERRQ(ierr);
81 
82   ierr = MatMarkDiagonal_SeqBAIJ(mpimat->A);CHKERRQ(ierr);
83   for(i=0;i<lm/bs;i++){
84     d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
85     o_nnz[i] = Ba->i[i+1] - Ba->i[i];
86   }
87 
88   ierr = MatCreate(((PetscObject)A)->comm,&M);CHKERRQ(ierr);
89   ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
90   ierr = MatSetType(M,newtype);CHKERRQ(ierr);
91   ierr = MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
92 
93   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
94 
95   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
96   ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);
97   for(i=rstart;i<rend;i++){
98     ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
99     j = 0;
100     ierr = MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);CHKERRQ(ierr);
101     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
102   }
103   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
104   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105 
106   if (reuse == MAT_REUSE_MATRIX) {
107     ierr = MatHeaderReplace(A,M);CHKERRQ(ierr);
108   } else {
109     *newmat = M;
110   }
111   PetscFunctionReturn(0);
112 }
113 EXTERN_C_END
114