xref: /petsc/src/mat/impls/shell/shellcnv.c (revision efe6450ae345901f2e8b07d9b8db05e083e9ced5)
173f4d377SMatthew Knepley /*$Id: shellcnv.c,v 1.20 2001/08/07 03:02:51 balay Exp $*/
26098dc10SBarry Smith 
3e090d566SSatish Balay #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
46098dc10SBarry Smith #include "src/vec/vecimpl.h"
56098dc10SBarry Smith 
64a2ae208SSatish Balay #undef __FUNCT__
74a2ae208SSatish Balay #define __FUNCT__ "MatConvert_Shell"
8676c34cdSKris Buschelman int MatConvert_Shell(Mat oldmat,MatType newtype,Mat *newmat) {
9676c34cdSKris Buschelman   Mat           mat;
106098dc10SBarry Smith   Vec           in,out;
11273d9f13SBarry Smith   int           ierr,i,M,m,*rows,start,end;
126098dc10SBarry Smith   MPI_Comm      comm;
1387828ca2SBarry Smith   PetscScalar   *array,zero = 0.0,one = 1.0;
146098dc10SBarry Smith 
153a40ed3dSBarry Smith   PetscFunctionBegin;
166098dc10SBarry Smith   comm = oldmat->comm;
176098dc10SBarry Smith 
186098dc10SBarry Smith   ierr = MatGetOwnershipRange(oldmat,&start,&end);CHKERRQ(ierr);
196098dc10SBarry Smith   ierr = VecCreateMPI(comm,end-start,PETSC_DECIDE,&in);CHKERRQ(ierr);
206098dc10SBarry Smith   ierr = VecDuplicate(in,&out);CHKERRQ(ierr);
216098dc10SBarry Smith   ierr = VecGetSize(in,&M);CHKERRQ(ierr);
226098dc10SBarry Smith   ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr);
23b0a32e0cSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(int),&rows);CHKERRQ(ierr);
246098dc10SBarry Smith   for (i=0; i<m; i++) {rows[i] = start + i;}
256098dc10SBarry Smith 
26676c34cdSKris Buschelman   ierr = MatCreate(comm,m,M,M,M,&mat);CHKERRQ(ierr);
27676c34cdSKris Buschelman   ierr = MatSetType(mat,newtype);CHKERRQ(ierr);
286098dc10SBarry Smith 
296098dc10SBarry Smith   for (i=0; i<M; i++) {
306098dc10SBarry Smith     ierr = VecSet(&zero,in);CHKERRQ(ierr);
316098dc10SBarry Smith     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
326098dc10SBarry Smith     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
336098dc10SBarry Smith     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
346098dc10SBarry Smith 
356098dc10SBarry Smith     ierr = MatMult(oldmat,in,out);CHKERRQ(ierr);
366098dc10SBarry Smith 
376098dc10SBarry Smith     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
38676c34cdSKris Buschelman     ierr = MatSetValues(mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
396098dc10SBarry Smith     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
406098dc10SBarry Smith 
416098dc10SBarry Smith   }
42606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
436098dc10SBarry Smith   ierr = VecDestroy(in);CHKERRQ(ierr);
446098dc10SBarry Smith   ierr = VecDestroy(out);CHKERRQ(ierr);
45676c34cdSKris Buschelman   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
46676c34cdSKris Buschelman   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
47676c34cdSKris Buschelman 
48676c34cdSKris Buschelman   /* Fake support for "inplace" convert. */
49*efe6450aSKris Buschelman   if (*newmat == oldmat) {
50*efe6450aSKris Buschelman     ierr = MatDestroy(oldmat);CHKERRQ(ierr);
51676c34cdSKris Buschelman   }
52676c34cdSKris Buschelman   *newmat = mat;
533a40ed3dSBarry Smith   PetscFunctionReturn(0);
546098dc10SBarry Smith }
556098dc10SBarry Smith 
566098dc10SBarry Smith 
576098dc10SBarry Smith 
58