xref: /petsc/src/mat/impls/sbaij/mpi/mpiaijsbaij.c (revision 1667e6754bbf5beed23274f5aabc810cb13136ba)
1 #define PETSCMAT_DLL
2 
3 #include "mpisbaij.h" /*I "petscmat.h" I*/
4 #include "../src/mat/impls/aij/mpi/mpiaij.h"
5 #include "private/matimpl.h"
6 #include "petscmat.h"
7 
8 EXTERN_C_BEGIN
9 #undef __FUNCT__
10 #define __FUNCT__ "MatConvert_MPIAIJ_MPISBAIJ"
11 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MPIAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
12 {
13   PetscErrorCode     ierr;
14   Mat                M;
15   PetscInt           *d_nnz,*o_nnz;
16   PetscInt           i,j,k,nz;
17   PetscInt           m,n,lm,ln;
18   PetscInt           rstart,rend;
19   const PetscScalar  *vwork;
20   const PetscInt     *cwork;
21 
22   PetscFunctionBegin;
23   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
24   ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr);
25   ierr = PetscMalloc2(lm*sizeof(PetscInt),PetscInt,&d_nnz,lm*sizeof(PetscInt),PetscInt,&o_nnz);CHKERRQ(ierr);
26 
27   k = 0;
28   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
29   for(i=rstart;i<rend;i++){
30     ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
31     d_nnz[k] = o_nnz[k] = 0;
32     for(j=0;j<nz;j++){
33       if(cwork[j] >= rstart && cwork[j] <= rend-1){ /* diagonal portion */
34 	if(cwork[j] > i)
35 	  d_nnz[k] += 1;
36       }
37       else{
38 	/* insert values only in upper triangular portion */
39 	if(cwork[j] >= rend)
40 	    o_nnz[k] += 1;  /* off-diagonal portion */
41       }
42     }
43     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
44     k++;
45   }
46 
47   ierr = MatCreate(((PetscObject)A)->comm,&M);CHKERRQ(ierr);
48   ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
49   ierr = MatSetType(M,newtype);CHKERRQ(ierr);
50   ierr = MatMPISBAIJSetPreallocation(M,1,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
51 
52   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
53 
54   for(i=rstart;i<rend;i++){
55     ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
56     for(j=0;j<nz;j++){
57       if(cwork[j] >= rstart && cwork[j] <= rend-1){ /* diagonal portion */
58 	if(cwork[j] > i)
59 	  ierr = MatSetValues(M,1,&i,1,&cwork[j],&vwork[j],INSERT_VALUES);CHKERRQ(ierr);
60       }else{
61 	/* insert values only in upper triangular portion */
62 	if(cwork[j] >= rend){
63 	  ierr = MatSetValues(M,1,&i,1,&cwork[j],&vwork[j],INSERT_VALUES);CHKERRQ(ierr);
64 	}
65       }
66     }
67     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
68     k++;
69   }
70   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
72 
73   if (A->hermitian){
74     ierr = MatSetOption(M,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr);
75   }
76 
77   if (reuse == MAT_REUSE_MATRIX) {
78     ierr = MatHeaderReplace(A,M);CHKERRQ(ierr);
79   } else
80     *newmat = M;
81 
82   PetscFunctionReturn(0);
83 }
84 EXTERN_C_END
85