xref: /petsc/src/mat/impls/shell/shellcnv.c (revision a7cc72af8dae79365b4ec1ba6d96bf0305bf6125)
1 
2 #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
3 #include "vecimpl.h"
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "MatConvert_Shell"
7 PetscErrorCode MatConvert_Shell(Mat oldmat,const MatType newtype,Mat *newmat)
8 {
9   Mat            mat;
10   Vec            in,out;
11   PetscErrorCode ierr;
12   PetscInt       i,M,m,*rows,start,end;
13   MPI_Comm       comm;
14   PetscScalar    *array,zero = 0.0,one = 1.0;
15 
16   PetscFunctionBegin;
17   comm = oldmat->comm;
18 
19   ierr = MatGetOwnershipRange(oldmat,&start,&end);CHKERRQ(ierr);
20   ierr = VecCreateMPI(comm,end-start,PETSC_DECIDE,&in);CHKERRQ(ierr);
21   ierr = VecDuplicate(in,&out);CHKERRQ(ierr);
22   ierr = VecGetSize(in,&M);CHKERRQ(ierr);
23   ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr);
24   ierr = PetscMalloc((m+1)*sizeof(int),&rows);CHKERRQ(ierr);
25   for (i=0; i<m; i++) {rows[i] = start + i;}
26 
27   ierr = MatCreate(comm,m,M,M,M,&mat);CHKERRQ(ierr);
28   ierr = MatSetType(mat,newtype);CHKERRQ(ierr);
29   ierr = MatSetBlockSize(mat,oldmat->bs);CHKERRQ(ierr);
30 
31   for (i=0; i<M; i++) {
32     ierr = VecSet(&zero,in);CHKERRQ(ierr);
33     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
34     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
35     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
36 
37     ierr = MatMult(oldmat,in,out);CHKERRQ(ierr);
38 
39     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
40     ierr = MatSetValues(mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
41     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
42 
43   }
44   ierr = PetscFree(rows);CHKERRQ(ierr);
45   ierr = VecDestroy(in);CHKERRQ(ierr);
46   ierr = VecDestroy(out);CHKERRQ(ierr);
47   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
49 
50   /* Fake support for "inplace" convert. */
51   if (*newmat == oldmat) {
52     ierr = MatDestroy(oldmat);CHKERRQ(ierr);
53   }
54   *newmat = mat;
55   PetscFunctionReturn(0);
56 }
57 
58 
59 
60