xref: /petsc/src/mat/impls/sbaij/mpi/mpiaijsbaij.c (revision 3c48a1e8da19189ff2402a4e41a2fc082d52c349)
1 
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"
6 
7 EXTERN_C_BEGIN
8 #undef __FUNCT__
9 #define __FUNCT__ "MatConvert_MPIAIJ_MPISBAIJ"
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;
22 
23   PetscFunctionBegin;
24   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
25   ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr);
26   ierr = PetscMalloc2(lm*sizeof(PetscInt),PetscInt,&d_nnz,lm*sizeof(PetscInt),PetscInt,&o_nnz);CHKERRQ(ierr);
27 
28   for(i=0;i<lm;i++){
29     d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
30     o_nnz[i] = Ba->i[i+1] - Ba->i[i];
31   }
32 
33   ierr = MatCreate(((PetscObject)A)->comm,&M);CHKERRQ(ierr);
34   ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
35   ierr = MatSetType(M,newtype);CHKERRQ(ierr);
36   ierr = MatMPISBAIJSetPreallocation(M,1,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
37 
38   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
39 
40   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
41   for(i=rstart;i<rend;i++){
42     ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
43     j = 0;
44     while (cwork[j] < i){ j++; nz--;}
45     ierr = MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);CHKERRQ(ierr);
46     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
47   }
48   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
49   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
50 
51   if (A->hermitian){
52     ierr = MatSetOption(M,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr);
53   }
54 
55   if (reuse == MAT_REUSE_MATRIX) {
56     ierr = MatHeaderReplace(A,M);CHKERRQ(ierr);
57   } else {
58     *newmat = M;
59   }
60   PetscFunctionReturn(0);
61 }
62 EXTERN_C_END
63