Actual source code: mpiaijsbaij.c
2: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> /*I "petscmat.h" I*/
3: #include <../src/mat/impls/aij/mpi/mpiaij.h>
4: #include <private/matimpl.h>
5: #include <petscmat.h>
10: PetscErrorCode MatConvert_MPIAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
11: {
12: PetscErrorCode ierr;
13: Mat M;
14: Mat_MPIAIJ *mpimat = (Mat_MPIAIJ*)A->data;
15: Mat_SeqAIJ *Aa = (Mat_SeqAIJ*)mpimat->A->data,*Ba = (Mat_SeqAIJ*)mpimat->B->data;
16: PetscInt *d_nnz,*o_nnz;
17: PetscInt i,j,nz;
18: PetscInt m,n,lm,ln;
19: PetscInt rstart,rend;
20: const PetscScalar *vwork;
21: const PetscInt *cwork;
24: if (!A->symmetric) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
25: MatGetSize(A,&m,&n);
26: MatGetLocalSize(A,&lm,&ln);
27: PetscMalloc2(lm,PetscInt,&d_nnz,lm,PetscInt,&o_nnz);
29: MatMarkDiagonal_SeqAIJ(mpimat->A);
30: for(i=0;i<lm;i++){
31: d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
32: o_nnz[i] = Ba->i[i+1] - Ba->i[i];
33: }
35: MatCreate(((PetscObject)A)->comm,&M);
36: MatSetSizes(M,lm,ln,m,n);
37: MatSetType(M,MATMPISBAIJ);
38: MatSeqSBAIJSetPreallocation(M,1,0,d_nnz);
39: MatMPISBAIJSetPreallocation(M,1,0,d_nnz,0,o_nnz);
41: PetscFree2(d_nnz,o_nnz);
43: MatGetOwnershipRange(A,&rstart,&rend);
44: for(i=rstart;i<rend;i++){
45: MatGetRow(A,i,&nz,&cwork,&vwork);
46: j = 0;
47: while (cwork[j] < i){ j++; nz--;}
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_REUSE_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 */
64: PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
65: {
66: PetscErrorCode ierr;
67: Mat M;
68: Mat_MPIBAIJ *mpimat = (Mat_MPIBAIJ*)A->data;
69: Mat_SeqBAIJ *Aa = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data;
70: PetscInt *d_nnz,*o_nnz;
71: PetscInt i,j,nz;
72: PetscInt m,n,lm,ln;
73: PetscInt rstart,rend;
74: const PetscScalar *vwork;
75: const PetscInt *cwork;
76: PetscInt bs = A->rmap->bs;
79: MatGetSize(A,&m,&n);
80: MatGetLocalSize(A,&lm,&ln);
81: PetscMalloc2(lm/bs,PetscInt,&d_nnz,lm/bs,PetscInt,&o_nnz);
82:
83: MatMarkDiagonal_SeqBAIJ(mpimat->A);
84: for(i=0;i<lm/bs;i++){
85: d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
86: o_nnz[i] = Ba->i[i+1] - Ba->i[i];
87: }
89: MatCreate(((PetscObject)A)->comm,&M);
90: MatSetSizes(M,lm,ln,m,n);
91: MatSetType(M,MATMPISBAIJ);
92: MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);
93: MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);
95: PetscFree2(d_nnz,o_nnz);
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_REUSE_MATRIX) {
109: MatHeaderReplace(A,M);
110: } else {
111: *newmat = M;
112: }
113: return(0);
114: }