xref: /petsc/src/mat/impls/shell/shellcnv.c (revision 3c94ec114b93ba0cfbf7331735f0e93b9db506db)
1 /*$Id: shellcnv.c,v 1.20 2001/08/07 03:02:51 balay Exp $*/
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 int MatConvert_Shell(Mat oldmat,const MatType newtype,Mat *newmat) {
9   Mat           mat;
10   Vec           in,out;
11   int           ierr,i,M,m,*rows,start,end;
12   MPI_Comm      comm;
13   PetscScalar   *array,zero = 0.0,one = 1.0;
14 
15   PetscFunctionBegin;
16   comm = oldmat->comm;
17 
18   ierr = MatGetOwnershipRange(oldmat,&start,&end);CHKERRQ(ierr);
19   ierr = VecCreateMPI(comm,end-start,PETSC_DECIDE,&in);CHKERRQ(ierr);
20   ierr = VecDuplicate(in,&out);CHKERRQ(ierr);
21   ierr = VecGetSize(in,&M);CHKERRQ(ierr);
22   ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr);
23   ierr = PetscMalloc((m+1)*sizeof(int),&rows);CHKERRQ(ierr);
24   for (i=0; i<m; i++) {rows[i] = start + i;}
25 
26   ierr = MatCreate(comm,m,M,M,M,&mat);CHKERRQ(ierr);
27   ierr = MatSetType(mat,newtype);CHKERRQ(ierr);
28 
29   for (i=0; i<M; i++) {
30     ierr = VecSet(&zero,in);CHKERRQ(ierr);
31     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
32     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
33     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
34 
35     ierr = MatMult(oldmat,in,out);CHKERRQ(ierr);
36 
37     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
38     ierr = MatSetValues(mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
39     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
40 
41   }
42   ierr = PetscFree(rows);CHKERRQ(ierr);
43   ierr = VecDestroy(in);CHKERRQ(ierr);
44   ierr = VecDestroy(out);CHKERRQ(ierr);
45   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
46   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
47 
48   /* Fake support for "inplace" convert. */
49   if (*newmat == oldmat) {
50     ierr = MatDestroy(oldmat);CHKERRQ(ierr);
51   }
52   *newmat = mat;
53   PetscFunctionReturn(0);
54 }
55 
56 
57 
58