Actual source code: mpiaijsbaij.c
petsc-3.9.4 2018-09-11
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: MatGetSize(A,&m,&n);
23: MatGetLocalSize(A,&lm,&ln);
24: PetscMalloc2(lm,&d_nnz,lm,&o_nnz);
26: MatMarkDiagonal_SeqAIJ(mpimat->A);
27: for (i=0; i<lm; i++) {
28: d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
29: o_nnz[i] = Ba->i[i+1] - Ba->i[i];
30: }
32: MatCreate(PetscObjectComm((PetscObject)A),&M);
33: MatSetSizes(M,lm,ln,m,n);
34: MatSetType(M,MATMPISBAIJ);
35: MatSeqSBAIJSetPreallocation(M,1,0,d_nnz);
36: MatMPISBAIJSetPreallocation(M,1,0,d_nnz,0,o_nnz);
38: PetscFree2(d_nnz,o_nnz);
40: MatGetOwnershipRange(A,&rstart,&rend);
41: for (i=rstart; i<rend; i++) {
42: MatGetRow(A,i,&nz,&cwork,&vwork);
43: j = 0;
44: while (cwork[j] < i) {
45: j++; nz--;
46: }
47: MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);
48: MatRestoreRow(A,i,&nz,&cwork,&vwork);
49: }
50: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
51: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
53: if (reuse == MAT_INPLACE_MATRIX) {
54: MatHeaderReplace(A,&M);
55: } else {
56: *newmat = M;
57: }
58: return(0);
59: }
60: /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */
61: PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
62: {
63: PetscErrorCode ierr;
64: Mat M;
65: Mat_MPIBAIJ *mpimat = (Mat_MPIBAIJ*)A->data;
66: Mat_SeqBAIJ *Aa = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data;
67: PetscInt *d_nnz,*o_nnz;
68: PetscInt i,j,nz;
69: PetscInt m,n,lm,ln;
70: PetscInt rstart,rend;
71: const PetscScalar *vwork;
72: const PetscInt *cwork;
73: PetscInt bs = A->rmap->bs;
76: MatGetSize(A,&m,&n);
77: MatGetLocalSize(A,&lm,&ln);
78: PetscMalloc2(lm/bs,&d_nnz,lm/bs,&o_nnz);
80: MatMarkDiagonal_SeqBAIJ(mpimat->A);
81: for (i=0; i<lm/bs; i++) {
82: d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
83: o_nnz[i] = Ba->i[i+1] - Ba->i[i];
84: }
86: MatCreate(PetscObjectComm((PetscObject)A),&M);
87: MatSetSizes(M,lm,ln,m,n);
88: MatSetType(M,MATMPISBAIJ);
89: MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);
90: MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);
92: PetscFree2(d_nnz,o_nnz);
94: MatGetOwnershipRange(A,&rstart,&rend);
95: MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
96: for (i=rstart; i<rend; i++) {
97: MatGetRow(A,i,&nz,&cwork,&vwork);
98: j = 0;
99: MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);
100: MatRestoreRow(A,i,&nz,&cwork,&vwork);
101: }
102: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
103: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
105: if (reuse == MAT_INPLACE_MATRIX) {
106: MatHeaderReplace(A,&M);
107: } else {
108: *newmat = M;
109: }
110: return(0);
111: }