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