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