Actual source code: convert.c
2: #include <private/matimpl.h>
6: /*
7: MatConvert_Basic - Converts from any input format to another format. For
8: parallel formats, the new matrix distribution is determined by PETSc.
10: Does not do preallocation so in general will be slow
11: */
12: PetscErrorCode MatConvert_Basic(Mat mat, const MatType newtype,MatReuse reuse,Mat *newmat)
13: {
14: Mat M;
15: const PetscScalar *vwork;
16: PetscErrorCode ierr;
17: PetscInt i,nz,m,n,rstart,rend,lm,ln;
18: const PetscInt *cwork;
19: PetscBool isseqsbaij,ismpisbaij;
22: MatGetSize(mat,&m,&n);
23: MatGetLocalSize(mat,&lm,&ln);
25: if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */
26:
27: MatCreate(((PetscObject)mat)->comm,&M);
28: MatSetSizes(M,lm,ln,m,n);
29: MatSetType(M,newtype);
30: PetscTypeCompare((PetscObject)M,MATSEQSBAIJ,&isseqsbaij);
31: PetscTypeCompare((PetscObject)M,MATMPISBAIJ,&ismpisbaij);
32: if (isseqsbaij || ismpisbaij) {MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);}
34: MatGetOwnershipRange(mat,&rstart,&rend);
35: for (i=rstart; i<rend; i++) {
36: MatGetRow(mat,i,&nz,&cwork,&vwork);
37: MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);
38: MatRestoreRow(mat,i,&nz,&cwork,&vwork);
39: }
40: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
41: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
43: if (reuse == MAT_REUSE_MATRIX) {
44: MatHeaderReplace(mat,M);
45: } else {
46: *newmat = M;
47: }
48: return(0);
49: }