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