xref: /petsc/src/mat/impls/sbaij/mpi/mpiaijsbaij.c (revision b2ccae6bdc8edea944f1c160ca3b2eb32c69ecb2)
1 #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> /*I "petscmat.h" I*/
2 #include <../src/mat/impls/aij/mpi/mpiaij.h>
3 #include <petsc/private/matimpl.h>
4 #include <petscmat.h>
5 
6 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ_Preallocate(Mat, PetscInt **);
7 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ_Preallocate(Mat, PetscInt **);
8 
9 /* The code is virtually identical to MatConvert_MPIAIJ_MPIBAIJ() */
10 PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_MPISBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
11 {
12   Mat         M;
13   Mat_MPIAIJ *mpimat = (Mat_MPIAIJ *)A->data;
14   PetscInt   *d_nnz, *o_nnz;
15   PetscInt    m, n, lm, ln, bs = A->rmap->bs;
16 
17   PetscFunctionBegin;
18   if (reuse != MAT_REUSE_MATRIX) {
19     PetscCall(MatDisAssemble_MPIAIJ(A, PETSC_FALSE));
20     PetscCall(MatGetSize(A, &m, &n));
21     PetscCall(MatGetLocalSize(A, &lm, &ln));
22     PetscCall(MatConvert_SeqAIJ_SeqSBAIJ_Preallocate(mpimat->A, &d_nnz));
23     PetscCall(MatConvert_SeqAIJ_SeqBAIJ_Preallocate(mpimat->B, &o_nnz));
24     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &M));
25     PetscCall(MatSetSizes(M, lm, ln, m, n));
26     PetscCall(MatSetType(M, MATMPISBAIJ));
27     PetscCall(MatSeqSBAIJSetPreallocation(M, bs, 0, d_nnz));
28     PetscCall(MatMPISBAIJSetPreallocation(M, bs, 0, d_nnz, 0, o_nnz));
29     PetscCall(PetscFree(d_nnz));
30     PetscCall(PetscFree(o_nnz));
31     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
32     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
33   } else M = *newmat;
34 
35   /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */
36   /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before                                          */
37   /*                      MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below                                        */
38   PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &M));
39 
40   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &M));
41   else *newmat = M;
42   PetscFunctionReturn(PETSC_SUCCESS);
43 }
44 
45 /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */
46 PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
47 {
48   Mat                M;
49   Mat_MPIBAIJ       *mpimat = (Mat_MPIBAIJ *)A->data;
50   Mat_SeqBAIJ       *Aa = (Mat_SeqBAIJ *)mpimat->A->data, *Ba = (Mat_SeqBAIJ *)mpimat->B->data;
51   PetscInt          *d_nnz, *o_nnz;
52   PetscInt           i, nz;
53   PetscInt           m, n, lm, ln;
54   PetscInt           rstart, rend;
55   const PetscScalar *vwork;
56   const PetscInt    *cwork;
57   PetscInt           bs = A->rmap->bs;
58 
59   PetscFunctionBegin;
60   if (reuse != MAT_REUSE_MATRIX) {
61     PetscCall(MatGetSize(A, &m, &n));
62     PetscCall(MatGetLocalSize(A, &lm, &ln));
63     PetscCall(PetscMalloc2(lm / bs, &d_nnz, lm / bs, &o_nnz));
64 
65     PetscCall(MatMarkDiagonal_SeqBAIJ(mpimat->A));
66     for (i = 0; i < lm / bs; i++) {
67       d_nnz[i] = Aa->i[i + 1] - Aa->diag[i];
68       o_nnz[i] = Ba->i[i + 1] - Ba->i[i];
69     }
70 
71     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &M));
72     PetscCall(MatSetSizes(M, lm, ln, m, n));
73     PetscCall(MatSetType(M, MATMPISBAIJ));
74     PetscCall(MatSeqSBAIJSetPreallocation(M, bs, 0, d_nnz));
75     PetscCall(MatMPISBAIJSetPreallocation(M, bs, 0, d_nnz, 0, o_nnz));
76 
77     PetscCall(PetscFree2(d_nnz, o_nnz));
78   } else M = *newmat;
79 
80   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
81   PetscCall(MatSetOption(M, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
82   for (i = rstart; i < rend; i++) {
83     PetscCall(MatGetRow(A, i, &nz, &cwork, &vwork));
84     PetscCall(MatSetValues(M, 1, &i, nz, cwork, vwork, INSERT_VALUES));
85     PetscCall(MatRestoreRow(A, i, &nz, &cwork, &vwork));
86   }
87   PetscCall(MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY));
88   PetscCall(MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY));
89 
90   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &M));
91   else *newmat = M;
92   PetscFunctionReturn(PETSC_SUCCESS);
93 }
94