xref: /petsc/src/mat/impls/sbaij/mpi/mpiaijsbaij.c (revision 1873fc9fedfc3bd8fc2f82ccd29ab8b44a80f1db)
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,bs=PetscAbs(A->rmap->bs);
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/bs,&d_nnz,lm/bs,&o_nnz);CHKERRQ(ierr);
26 
27     ierr = MatMarkDiagonal_SeqAIJ(mpimat->A);CHKERRQ(ierr);
28     for (i=0; i<lm/bs; i++) {
29       d_nnz[i] = (Aa->i[i*bs+1] - Aa->diag[i*bs])/bs;
30       o_nnz[i] = (Ba->i[i*bs+1] - Ba->i[i*bs])/bs;
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,bs,0,d_nnz);CHKERRQ(ierr);
37     ierr = MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
38     ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
39   } else M = *newmat;
40 
41   if (bs == 1) {
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   } else {
55     ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);
56     /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */
57     /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before                                          */
58     /*                      MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below                                        */
59     ierr = MatConvert_Basic(A,newtype,MAT_REUSE_MATRIX,&M);CHKERRQ(ierr);
60   }
61 
62   if (reuse == MAT_INPLACE_MATRIX) {
63     ierr = MatHeaderReplace(A,&M);CHKERRQ(ierr);
64   } else *newmat = M;
65   PetscFunctionReturn(0);
66 }
67 /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */
68 PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
69 {
70   PetscErrorCode    ierr;
71   Mat               M;
72   Mat_MPIBAIJ       *mpimat = (Mat_MPIBAIJ*)A->data;
73   Mat_SeqBAIJ       *Aa     = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data;
74   PetscInt          *d_nnz,*o_nnz;
75   PetscInt          i,nz;
76   PetscInt          m,n,lm,ln;
77   PetscInt          rstart,rend;
78   const PetscScalar *vwork;
79   const PetscInt    *cwork;
80   PetscInt          bs = A->rmap->bs;
81 
82   PetscFunctionBegin;
83   if (reuse != MAT_REUSE_MATRIX) {
84     ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
85     ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr);
86     ierr = PetscMalloc2(lm/bs,&d_nnz,lm/bs,&o_nnz);CHKERRQ(ierr);
87 
88     ierr = MatMarkDiagonal_SeqBAIJ(mpimat->A);CHKERRQ(ierr);
89     for (i=0; i<lm/bs; i++) {
90       d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
91       o_nnz[i] = Ba->i[i+1] - Ba->i[i];
92     }
93 
94     ierr = MatCreate(PetscObjectComm((PetscObject)A),&M);CHKERRQ(ierr);
95     ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
96     ierr = MatSetType(M,MATMPISBAIJ);CHKERRQ(ierr);
97     ierr = MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);CHKERRQ(ierr);
98     ierr = MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
99 
100     ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
101   } else M = *newmat;
102 
103   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
104   ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);
105   for (i=rstart; i<rend; i++) {
106     ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
107     ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
108     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
109   }
110   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
111   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
112 
113   if (reuse == MAT_INPLACE_MATRIX) {
114     ierr = MatHeaderReplace(A,&M);CHKERRQ(ierr);
115   } else *newmat = M;
116   PetscFunctionReturn(0);
117 }
118