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