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