Actual source code: mpiaijsbaij.c
petsc-3.14.6 2021-03-30
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 && !A->hermitian) SETERRQ(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)");
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: if (Aa->i[i*bs+1] == Aa->diag[i*bs]) { /* misses diagonal entry */
30: d_nnz[i] = (Aa->i[i*bs+1] - Aa->i[i*bs])/bs;
31: } else {
32: d_nnz[i] = (Aa->i[i*bs+1] - Aa->diag[i*bs])/bs;
33: }
34: o_nnz[i] = (Ba->i[i*bs+1] - Ba->i[i*bs])/bs;
35: }
37: MatCreate(PetscObjectComm((PetscObject)A),&M);
38: MatSetSizes(M,lm,ln,m,n);
39: MatSetType(M,MATMPISBAIJ);
40: MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);
41: MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);
42: PetscFree2(d_nnz,o_nnz);
43: } else M = *newmat;
45: if (bs == 1) {
46: MatGetOwnershipRange(A,&rstart,&rend);
47: for (i=rstart; i<rend; i++) {
48: MatGetRow(A,i,&nz,&cwork,&vwork);
49: if (nz) {
50: j = 0;
51: while (cwork[j] < i) {
52: j++; nz--;
53: }
54: MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);
55: }
56: MatRestoreRow(A,i,&nz,&cwork,&vwork);
57: }
58: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
59: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
60: } else {
61: 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: MatConvert_Basic(A,newtype,MAT_REUSE_MATRIX,&M);
66: }
68: if (reuse == MAT_INPLACE_MATRIX) {
69: MatHeaderReplace(A,&M);
70: } else *newmat = M;
71: return(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: PetscErrorCode ierr;
77: Mat M;
78: Mat_MPIBAIJ *mpimat = (Mat_MPIBAIJ*)A->data;
79: Mat_SeqBAIJ *Aa = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data;
80: PetscInt *d_nnz,*o_nnz;
81: PetscInt i,nz;
82: PetscInt m,n,lm,ln;
83: PetscInt rstart,rend;
84: const PetscScalar *vwork;
85: const PetscInt *cwork;
86: PetscInt bs = A->rmap->bs;
89: if (reuse != MAT_REUSE_MATRIX) {
90: MatGetSize(A,&m,&n);
91: MatGetLocalSize(A,&lm,&ln);
92: PetscMalloc2(lm/bs,&d_nnz,lm/bs,&o_nnz);
94: MatMarkDiagonal_SeqBAIJ(mpimat->A);
95: for (i=0; i<lm/bs; i++) {
96: d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
97: o_nnz[i] = Ba->i[i+1] - Ba->i[i];
98: }
100: MatCreate(PetscObjectComm((PetscObject)A),&M);
101: MatSetSizes(M,lm,ln,m,n);
102: MatSetType(M,MATMPISBAIJ);
103: MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);
104: MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);
106: PetscFree2(d_nnz,o_nnz);
107: } else M = *newmat;
109: MatGetOwnershipRange(A,&rstart,&rend);
110: MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
111: for (i=rstart; i<rend; i++) {
112: MatGetRow(A,i,&nz,&cwork,&vwork);
113: MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);
114: MatRestoreRow(A,i,&nz,&cwork,&vwork);
115: }
116: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
117: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
119: if (reuse == MAT_INPLACE_MATRIX) {
120: MatHeaderReplace(A,&M);
121: } else *newmat = M;
122: return(0);
123: }