xref: /petsc/src/mat/impls/sbaij/mpi/mpiaijsbaij.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
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   Mat                M;
9   Mat_MPIAIJ        *mpimat = (Mat_MPIAIJ *)A->data;
10   Mat_SeqAIJ        *Aa = (Mat_SeqAIJ *)mpimat->A->data, *Ba = (Mat_SeqAIJ *)mpimat->B->data;
11   PetscInt          *d_nnz, *o_nnz;
12   PetscInt           i, j, nz;
13   PetscInt           m, n, lm, ln;
14   PetscInt           rstart, rend, bs = PetscAbs(A->rmap->bs);
15   const PetscScalar *vwork;
16   const PetscInt    *cwork;
17 
18   PetscFunctionBegin;
19   PetscCheck(A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Matrix must be symmetric or Hermitian. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE) or MatSetOption(mat,MAT_HERMITIAN,PETSC_TRUE)");
20   if (reuse != MAT_REUSE_MATRIX) {
21     PetscCall(MatGetSize(A, &m, &n));
22     PetscCall(MatGetLocalSize(A, &lm, &ln));
23     PetscCall(PetscMalloc2(lm / bs, &d_nnz, lm / bs, &o_nnz));
24 
25     PetscCall(MatMarkDiagonal_SeqAIJ(mpimat->A));
26     for (i = 0; i < lm / bs; i++) {
27       if (Aa->i[i * bs + 1] == Aa->diag[i * bs]) { /* misses diagonal entry */
28         d_nnz[i] = (Aa->i[i * bs + 1] - Aa->i[i * bs]) / bs;
29       } else {
30         d_nnz[i] = (Aa->i[i * bs + 1] - Aa->diag[i * bs]) / bs;
31       }
32       o_nnz[i] = (Ba->i[i * bs + 1] - Ba->i[i * bs]) / bs;
33     }
34 
35     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &M));
36     PetscCall(MatSetSizes(M, lm, ln, m, n));
37     PetscCall(MatSetType(M, MATMPISBAIJ));
38     PetscCall(MatSeqSBAIJSetPreallocation(M, bs, 0, d_nnz));
39     PetscCall(MatMPISBAIJSetPreallocation(M, bs, 0, d_nnz, 0, o_nnz));
40     PetscCall(PetscFree2(d_nnz, o_nnz));
41   } else M = *newmat;
42 
43   if (bs == 1) {
44     PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
45     for (i = rstart; i < rend; i++) {
46       PetscCall(MatGetRow(A, i, &nz, &cwork, &vwork));
47       if (nz) {
48         j = 0;
49         while (cwork[j] < i) {
50           j++;
51           nz--;
52         }
53         PetscCall(MatSetValues(M, 1, &i, nz, cwork + j, vwork + j, INSERT_VALUES));
54       }
55       PetscCall(MatRestoreRow(A, i, &nz, &cwork, &vwork));
56     }
57     PetscCall(MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY));
58     PetscCall(MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY));
59   } else {
60     PetscCall(MatSetOption(M, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
61     /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */
62     /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before                                          */
63     /*                      MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below                                        */
64     PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &M));
65   }
66 
67   if (reuse == MAT_INPLACE_MATRIX) {
68     PetscCall(MatHeaderReplace(A, &M));
69   } else *newmat = M;
70   PetscFunctionReturn(0);
71 }
72 /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */
73 PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) {
74   Mat                M;
75   Mat_MPIBAIJ       *mpimat = (Mat_MPIBAIJ *)A->data;
76   Mat_SeqBAIJ       *Aa = (Mat_SeqBAIJ *)mpimat->A->data, *Ba = (Mat_SeqBAIJ *)mpimat->B->data;
77   PetscInt          *d_nnz, *o_nnz;
78   PetscInt           i, nz;
79   PetscInt           m, n, lm, ln;
80   PetscInt           rstart, rend;
81   const PetscScalar *vwork;
82   const PetscInt    *cwork;
83   PetscInt           bs = A->rmap->bs;
84 
85   PetscFunctionBegin;
86   if (reuse != MAT_REUSE_MATRIX) {
87     PetscCall(MatGetSize(A, &m, &n));
88     PetscCall(MatGetLocalSize(A, &lm, &ln));
89     PetscCall(PetscMalloc2(lm / bs, &d_nnz, lm / bs, &o_nnz));
90 
91     PetscCall(MatMarkDiagonal_SeqBAIJ(mpimat->A));
92     for (i = 0; i < lm / bs; i++) {
93       d_nnz[i] = Aa->i[i + 1] - Aa->diag[i];
94       o_nnz[i] = Ba->i[i + 1] - Ba->i[i];
95     }
96 
97     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &M));
98     PetscCall(MatSetSizes(M, lm, ln, m, n));
99     PetscCall(MatSetType(M, MATMPISBAIJ));
100     PetscCall(MatSeqSBAIJSetPreallocation(M, bs, 0, d_nnz));
101     PetscCall(MatMPISBAIJSetPreallocation(M, bs, 0, d_nnz, 0, o_nnz));
102 
103     PetscCall(PetscFree2(d_nnz, o_nnz));
104   } else M = *newmat;
105 
106   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
107   PetscCall(MatSetOption(M, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
108   for (i = rstart; i < rend; i++) {
109     PetscCall(MatGetRow(A, i, &nz, &cwork, &vwork));
110     PetscCall(MatSetValues(M, 1, &i, nz, cwork, vwork, INSERT_VALUES));
111     PetscCall(MatRestoreRow(A, i, &nz, &cwork, &vwork));
112   }
113   PetscCall(MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY));
114   PetscCall(MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY));
115 
116   if (reuse == MAT_INPLACE_MATRIX) {
117     PetscCall(MatHeaderReplace(A, &M));
118   } else *newmat = M;
119   PetscFunctionReturn(0);
120 }
121