xref: /petsc/src/mat/utils/convert.c (revision 1a0cae60fb91c8e2ce9276d969629479e2ae5b13)
1 
2 #include <petsc-private/matimpl.h>
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "MatConvert_Basic"
6 /*
7   MatConvert_Basic - Converts from any input format to another format. For
8   parallel formats, the new matrix distribution is determined by PETSc.
9 
10   Does not do preallocation so in general will be slow
11  */
12 PetscErrorCode MatConvert_Basic(Mat mat, MatType newtype,MatReuse reuse,Mat *newmat)
13 {
14   Mat               M;
15   const PetscScalar *vwork;
16   PetscErrorCode    ierr;
17   PetscInt          i,j,nz,m,n,rstart,rend,lm,ln,prbs,pcbs,cstart,cend,*dnz,*onz;
18   const PetscInt    *cwork;
19   PetscBool         isseqsbaij,ismpisbaij,isseqbaij,ismpibaij;
20 
21   PetscFunctionBegin;
22   ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr);
23   ierr = MatGetLocalSize(mat,&lm,&ln);CHKERRQ(ierr);
24 
25   if (ln == n) ln = PETSC_DECIDE; /* try to preserve column ownership */
26 
27   ierr = MatCreate(PetscObjectComm((PetscObject)mat),&M);CHKERRQ(ierr);
28   ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr);
29   ierr = MatSetBlockSizes(M,mat->rmap->bs,mat->cmap->bs);CHKERRQ(ierr);
30   ierr = MatSetType(M,newtype);CHKERRQ(ierr);
31   ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
32   ierr = PetscObjectTypeCompare((PetscObject)M,MATMPISBAIJ,&ismpisbaij);CHKERRQ(ierr);
33   if (isseqsbaij || ismpisbaij) {ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);}
34   ierr = PetscObjectTypeCompare((PetscObject)M,MATSEQBAIJ,&isseqbaij);CHKERRQ(ierr);
35   ierr = PetscObjectTypeCompare((PetscObject)M,MATMPIBAIJ,&ismpibaij);CHKERRQ(ierr);
36 
37   /* Preallocation block sizes.  (S)BAIJ matrices will have one index per block. */
38   prbs = (isseqbaij || ismpibaij || isseqsbaij || ismpisbaij) ? M->rmap->bs : 1;
39   pcbs = (isseqbaij || ismpibaij || isseqsbaij || ismpisbaij) ? M->cmap->bs : 1;
40 
41   ierr = PetscMalloc2(lm/prbs,PetscInt,&dnz,lm/prbs,PetscInt,&onz);CHKERRQ(ierr);
42   ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr);
43   ierr = MatGetOwnershipRangeColumn(mat,&cstart,&cend);CHKERRQ(ierr);
44   for (i=0; i<lm; i+=prbs) {
45     ierr = MatGetRow(mat,rstart+i,&nz,&cwork,NULL);CHKERRQ(ierr);
46     dnz[i] = 0;
47     onz[i] = 0;
48     for (j=0; j<nz; j+=pcbs) {
49       if ((isseqsbaij || ismpisbaij) && cwork[j] < rstart+i) continue;
50       if (cstart <= cwork[j] && cwork[j] < cend) dnz[i]++;
51       else                                       onz[i]++;
52     }
53     ierr = MatRestoreRow(mat,rstart+i,&nz,&cwork,NULL);CHKERRQ(ierr);
54   }
55   ierr = MatXAIJSetPreallocation(M,M->rmap->bs,dnz,onz,dnz,onz);CHKERRQ(ierr);
56   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
57 
58   for (i=rstart; i<rend; i++) {
59     ierr = MatGetRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
60     ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
61     ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
62   }
63   ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
64   ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
65 
66   if (reuse == MAT_REUSE_MATRIX) {
67     ierr = MatHeaderReplace(mat,M);CHKERRQ(ierr);
68   } else {
69     *newmat = M;
70   }
71   PetscFunctionReturn(0);
72 }
73