16098dc10SBarry Smith 2e090d566SSatish Balay #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 33c94ec11SBarry Smith #include "vecimpl.h" 46098dc10SBarry Smith 54a2ae208SSatish Balay #undef __FUNCT__ 64a2ae208SSatish Balay #define __FUNCT__ "MatConvert_Shell" 7ceb03754SKris Buschelman PetscErrorCode MatConvert_Shell(Mat oldmat,const MatType newtype,MatReuse reuse,Mat *newmat) 8dfbe8321SBarry Smith { 9676c34cdSKris Buschelman Mat mat; 106098dc10SBarry Smith Vec in,out; 11dfbe8321SBarry Smith PetscErrorCode ierr; 12b24ad042SBarry Smith PetscInt i,M,m,*rows,start,end; 136098dc10SBarry Smith MPI_Comm comm; 1487828ca2SBarry Smith PetscScalar *array,zero = 0.0,one = 1.0; 156098dc10SBarry Smith 163a40ed3dSBarry Smith PetscFunctionBegin; 176098dc10SBarry Smith comm = oldmat->comm; 186098dc10SBarry Smith 196098dc10SBarry Smith ierr = MatGetOwnershipRange(oldmat,&start,&end);CHKERRQ(ierr); 206098dc10SBarry Smith ierr = VecCreateMPI(comm,end-start,PETSC_DECIDE,&in);CHKERRQ(ierr); 216098dc10SBarry Smith ierr = VecDuplicate(in,&out);CHKERRQ(ierr); 226098dc10SBarry Smith ierr = VecGetSize(in,&M);CHKERRQ(ierr); 236098dc10SBarry Smith ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr); 24b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&rows);CHKERRQ(ierr); 256098dc10SBarry Smith for (i=0; i<m; i++) {rows[i] = start + i;} 266098dc10SBarry Smith 27676c34cdSKris Buschelman ierr = MatCreate(comm,m,M,M,M,&mat);CHKERRQ(ierr); 28676c34cdSKris Buschelman ierr = MatSetType(mat,newtype);CHKERRQ(ierr); 29521d7252SBarry Smith ierr = MatSetBlockSize(mat,oldmat->bs);CHKERRQ(ierr); 306098dc10SBarry Smith 316098dc10SBarry Smith for (i=0; i<M; i++) { 326098dc10SBarry Smith ierr = VecSet(&zero,in);CHKERRQ(ierr); 336098dc10SBarry Smith ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 346098dc10SBarry Smith ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 356098dc10SBarry Smith ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 366098dc10SBarry Smith 376098dc10SBarry Smith ierr = MatMult(oldmat,in,out);CHKERRQ(ierr); 386098dc10SBarry Smith 396098dc10SBarry Smith ierr = VecGetArray(out,&array);CHKERRQ(ierr); 40676c34cdSKris Buschelman ierr = MatSetValues(mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 416098dc10SBarry Smith ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 426098dc10SBarry Smith 436098dc10SBarry Smith } 44606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 456098dc10SBarry Smith ierr = VecDestroy(in);CHKERRQ(ierr); 466098dc10SBarry Smith ierr = VecDestroy(out);CHKERRQ(ierr); 47676c34cdSKris Buschelman ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 48676c34cdSKris Buschelman ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 49676c34cdSKris Buschelman 50ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 51*c3d102feSKris Buschelman ierr = MatHeaderCopy(oldmat,mat);CHKERRQ(ierr); 52*c3d102feSKris Buschelman } else { 53676c34cdSKris Buschelman *newmat = mat; 54*c3d102feSKris Buschelman } 553a40ed3dSBarry Smith PetscFunctionReturn(0); 566098dc10SBarry Smith } 576098dc10SBarry Smith 586098dc10SBarry Smith 596098dc10SBarry Smith 60