Actual source code: mpiaijsbaij.c
petsc-3.13.1 2020-05-02
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,bs=PetscAbs(A->rmap->bs);
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/bs,&d_nnz,lm/bs,&o_nnz);
27: MatMarkDiagonal_SeqAIJ(mpimat->A);
28: for (i=0; i<lm/bs; i++) {
29: d_nnz[i] = (Aa->i[i*bs+1] - Aa->diag[i*bs])/bs;
30: o_nnz[i] = (Ba->i[i*bs+1] - Ba->i[i*bs])/bs;
31: }
33: MatCreate(PetscObjectComm((PetscObject)A),&M);
34: MatSetSizes(M,lm,ln,m,n);
35: MatSetType(M,MATMPISBAIJ);
36: MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);
37: MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);
38: PetscFree2(d_nnz,o_nnz);
39: } else M = *newmat;
41: if (bs == 1) {
42: MatGetOwnershipRange(A,&rstart,&rend);
43: for (i=rstart; i<rend; i++) {
44: MatGetRow(A,i,&nz,&cwork,&vwork);
45: j = 0;
46: while (cwork[j] < i) {
47: j++; nz--;
48: }
49: MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);
50: MatRestoreRow(A,i,&nz,&cwork,&vwork);
51: }
52: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
53: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
54: } else {
55: MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
56: /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */
57: /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before */
58: /* MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below */
59: MatConvert_Basic(A,newtype,MAT_REUSE_MATRIX,&M);
60: }
62: if (reuse == MAT_INPLACE_MATRIX) {
63: MatHeaderReplace(A,&M);
64: } else *newmat = M;
65: return(0);
66: }
67: /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */
68: PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
69: {
70: PetscErrorCode ierr;
71: Mat M;
72: Mat_MPIBAIJ *mpimat = (Mat_MPIBAIJ*)A->data;
73: Mat_SeqBAIJ *Aa = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data;
74: PetscInt *d_nnz,*o_nnz;
75: PetscInt i,nz;
76: PetscInt m,n,lm,ln;
77: PetscInt rstart,rend;
78: const PetscScalar *vwork;
79: const PetscInt *cwork;
80: PetscInt bs = A->rmap->bs;
83: if (reuse != MAT_REUSE_MATRIX) {
84: MatGetSize(A,&m,&n);
85: MatGetLocalSize(A,&lm,&ln);
86: PetscMalloc2(lm/bs,&d_nnz,lm/bs,&o_nnz);
88: MatMarkDiagonal_SeqBAIJ(mpimat->A);
89: for (i=0; i<lm/bs; i++) {
90: d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
91: o_nnz[i] = Ba->i[i+1] - Ba->i[i];
92: }
94: MatCreate(PetscObjectComm((PetscObject)A),&M);
95: MatSetSizes(M,lm,ln,m,n);
96: MatSetType(M,MATMPISBAIJ);
97: MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);
98: MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);
100: PetscFree2(d_nnz,o_nnz);
101: } else M = *newmat;
103: MatGetOwnershipRange(A,&rstart,&rend);
104: MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
105: for (i=rstart; i<rend; i++) {
106: MatGetRow(A,i,&nz,&cwork,&vwork);
107: MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);
108: MatRestoreRow(A,i,&nz,&cwork,&vwork);
109: }
110: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
111: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
113: if (reuse == MAT_INPLACE_MATRIX) {
114: MatHeaderReplace(A,&M);
115: } else *newmat = M;
116: return(0);
117: }