xref: /petsc/src/mat/impls/sbaij/mpi/mpiaijsbaij.c (revision feff33ee0b5b037fa8f9f294dede656a2f85cc47)
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   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
23   ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr);
24   ierr = PetscMalloc2(lm,&d_nnz,lm,&o_nnz);CHKERRQ(ierr);
25 
26   ierr = MatMarkDiagonal_SeqAIJ(mpimat->A);CHKERRQ(ierr);
27   for (i=0; i<lm; i++) {
28     d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
29     o_nnz[i] = Ba->i[i+1] - Ba->i[i];
30   }
31 
32   ierr = MatCreate(PetscObjectComm((PetscObject)A),&M);CHKERRQ(ierr);
33   ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
34   ierr = MatSetType(M,MATMPISBAIJ);CHKERRQ(ierr);
35   ierr = MatSeqSBAIJSetPreallocation(M,1,0,d_nnz);CHKERRQ(ierr);
36   ierr = MatMPISBAIJSetPreallocation(M,1,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
37 
38   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
39 
40   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
41   for (i=rstart; i<rend; i++) {
42     ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
43     j    = 0;
44     while (cwork[j] < i) {
45       j++; nz--;
46     }
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_INPLACE_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 PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
62 {
63   PetscErrorCode    ierr;
64   Mat               M;
65   Mat_MPIBAIJ       *mpimat = (Mat_MPIBAIJ*)A->data;
66   Mat_SeqBAIJ       *Aa     = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data;
67   PetscInt          *d_nnz,*o_nnz;
68   PetscInt          i,j,nz;
69   PetscInt          m,n,lm,ln;
70   PetscInt          rstart,rend;
71   const PetscScalar *vwork;
72   const PetscInt    *cwork;
73   PetscInt          bs = A->rmap->bs;
74 
75   PetscFunctionBegin;
76   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
77   ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr);
78   ierr = PetscMalloc2(lm/bs,&d_nnz,lm/bs,&o_nnz);CHKERRQ(ierr);
79 
80   ierr = MatMarkDiagonal_SeqBAIJ(mpimat->A);CHKERRQ(ierr);
81   for (i=0; i<lm/bs; i++) {
82     d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
83     o_nnz[i] = Ba->i[i+1] - Ba->i[i];
84   }
85 
86   ierr = MatCreate(PetscObjectComm((PetscObject)A),&M);CHKERRQ(ierr);
87   ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
88   ierr = MatSetType(M,MATMPISBAIJ);CHKERRQ(ierr);
89   ierr = MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);CHKERRQ(ierr);
90   ierr = MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
91 
92   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
93 
94   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
95   ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);
96   for (i=rstart; i<rend; i++) {
97     ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
98     j    = 0;
99     ierr = MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);CHKERRQ(ierr);
100     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
101   }
102   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
103   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
104 
105   if (reuse == MAT_INPLACE_MATRIX) {
106     ierr = MatHeaderReplace(A,&M);CHKERRQ(ierr);
107   } else {
108     *newmat = M;
109   }
110   PetscFunctionReturn(0);
111 }
112