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