Actual source code: mpiaijsbaij.c
petsc-3.12.5 2020-03-29
2: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
3: #include <../src/mat/impls/aij/mpi/mpiaij.h>
4: #include <petsc/private/matimpl.h>
5: #include <petscmat.h>
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;
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: MatGetSize(A,&m,&n);
24: MatGetLocalSize(A,&lm,&ln);
25: PetscMalloc2(lm,&d_nnz,lm,&o_nnz);
27: MatMarkDiagonal_SeqAIJ(mpimat->A);
28: for (i=0; i<lm; i++) {
29: d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
30: o_nnz[i] = Ba->i[i+1] - Ba->i[i];
31: }
33: MatCreate(PetscObjectComm((PetscObject)A),&M);
34: MatSetSizes(M,lm,ln,m,n);
35: MatSetType(M,MATMPISBAIJ);
36: MatSeqSBAIJSetPreallocation(M,1,0,d_nnz);
37: MatMPISBAIJSetPreallocation(M,1,0,d_nnz,0,o_nnz);
38: PetscFree2(d_nnz,o_nnz);
39: } else M = *newmat;
41: MatGetOwnershipRange(A,&rstart,&rend);
42: for (i=rstart; i<rend; i++) {
43: MatGetRow(A,i,&nz,&cwork,&vwork);
44: j = 0;
45: while (cwork[j] < i) {
46: j++; nz--;
47: }
48: MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);
49: MatRestoreRow(A,i,&nz,&cwork,&vwork);
50: }
51: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
52: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
54: if (reuse == MAT_INPLACE_MATRIX) {
55: MatHeaderReplace(A,&M);
56: } else {
57: *newmat = M;
58: }
59: return(0);
60: }
61: /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */
62: PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
63: {
64: PetscErrorCode ierr;
65: Mat M;
66: Mat_MPIBAIJ *mpimat = (Mat_MPIBAIJ*)A->data;
67: Mat_SeqBAIJ *Aa = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data;
68: PetscInt *d_nnz,*o_nnz;
69: PetscInt i,j,nz;
70: PetscInt m,n,lm,ln;
71: PetscInt rstart,rend;
72: const PetscScalar *vwork;
73: const PetscInt *cwork;
74: PetscInt bs = A->rmap->bs;
77: if (reuse != MAT_REUSE_MATRIX) {
78: MatGetSize(A,&m,&n);
79: MatGetLocalSize(A,&lm,&ln);
80: PetscMalloc2(lm/bs,&d_nnz,lm/bs,&o_nnz);
82: MatMarkDiagonal_SeqBAIJ(mpimat->A);
83: for (i=0; i<lm/bs; i++) {
84: d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
85: o_nnz[i] = Ba->i[i+1] - Ba->i[i];
86: }
88: MatCreate(PetscObjectComm((PetscObject)A),&M);
89: MatSetSizes(M,lm,ln,m,n);
90: MatSetType(M,MATMPISBAIJ);
91: MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);
92: MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);
94: PetscFree2(d_nnz,o_nnz);
95: } else M = *newmat;
97: MatGetOwnershipRange(A,&rstart,&rend);
98: MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
99: for (i=rstart; i<rend; i++) {
100: MatGetRow(A,i,&nz,&cwork,&vwork);
101: j = 0;
102: MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);
103: MatRestoreRow(A,i,&nz,&cwork,&vwork);
104: }
105: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
106: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
108: if (reuse == MAT_INPLACE_MATRIX) {
109: MatHeaderReplace(A,&M);
110: } else {
111: *newmat = M;
112: }
113: return(0);
114: }