1*a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*a5eb4965SSatish Balay static char vcid[] = "$Id: convert.c,v 1.58 1997/01/27 18:17:25 bsmith Exp balay $"; 356fe5c5cSLois Curfman McInnes #endif 456fe5c5cSLois Curfman McInnes 594a9d846SBarry Smith #include "src/mat/matimpl.h" 68b6375c0SLois Curfman McInnes 75615d1e5SSatish Balay #undef __FUNC__ 85615d1e5SSatish Balay #define __FUNC__ "MatConvert_Basic" 98b6375c0SLois Curfman McInnes /* 108b6375c0SLois Curfman McInnes MatConvert_Basic - Converts from any input format to another format. For 118b6375c0SLois Curfman McInnes parallel formats, the new matrix distribution is determined by PETSc. 128b6375c0SLois Curfman McInnes */ 138b6375c0SLois Curfman McInnes int MatConvert_Basic(Mat mat,MatType newtype,Mat *M) 148b6375c0SLois Curfman McInnes { 158b6375c0SLois Curfman McInnes Scalar *vwork; 1625cdf11fSBarry Smith int ierr, i, nz, m, n, *cwork, rstart, rend,flg; 1725cdf11fSBarry Smith 188b6375c0SLois Curfman McInnes ierr = MatGetSize(mat,&m,&n); CHKERRQ(ierr); 197181dc05SLois Curfman McInnes if (newtype == MATSAME) newtype = (MatType)mat->type; 208b6375c0SLois Curfman McInnes switch (newtype) { 218b6375c0SLois Curfman McInnes case MATSEQAIJ: 22b4fd4287SBarry Smith ierr = MatCreateSeqAIJ(mat->comm,m,n,0,PETSC_NULL,M); CHKERRQ(ierr); 238b6375c0SLois Curfman McInnes break; 248b6375c0SLois Curfman McInnes case MATMPIROWBS: 25e3372554SBarry Smith if (m != n) SETERRQ(1,0,"MATMPIROWBS matrix must be square"); 263b2fbd54SBarry Smith ierr = MatCreateMPIRowbs(mat->comm,PETSC_DECIDE,m,0,PETSC_NULL, 27b4fd4287SBarry Smith PETSC_NULL,M); CHKERRQ(ierr); 288b6375c0SLois Curfman McInnes break; 298b6375c0SLois Curfman McInnes case MATMPIAIJ: 303b2fbd54SBarry Smith ierr = MatCreateMPIAIJ(mat->comm,PETSC_DECIDE,PETSC_DECIDE, 31b4fd4287SBarry Smith m,n,0,PETSC_NULL,0,PETSC_NULL,M); CHKERRQ(ierr); 328b6375c0SLois Curfman McInnes break; 338b6375c0SLois Curfman McInnes case MATSEQDENSE: 34b4fd4287SBarry Smith ierr = MatCreateSeqDense(mat->comm,m,n,PETSC_NULL,M); CHKERRQ(ierr); 358b6375c0SLois Curfman McInnes break; 368b6375c0SLois Curfman McInnes case MATMPIDENSE: 373b2fbd54SBarry Smith ierr = MatCreateMPIDense(mat->comm,PETSC_DECIDE,PETSC_DECIDE, 38b4fd4287SBarry Smith m,n,PETSC_NULL,M); CHKERRQ(ierr); 398b6375c0SLois Curfman McInnes break; 408b6375c0SLois Curfman McInnes case MATSEQBDIAG: 418b6375c0SLois Curfman McInnes { 42537820f0SBarry Smith int bs = 1; /* Default block size = 1 */ 43537820f0SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 44537820f0SBarry Smith ierr = MatCreateSeqBDiag(mat->comm,m,n,0,bs,PETSC_NULL,PETSC_NULL,M);CHKERRQ(ierr); 458b6375c0SLois Curfman McInnes break; 468b6375c0SLois Curfman McInnes } 478b6375c0SLois Curfman McInnes case MATMPIBDIAG: 488b6375c0SLois Curfman McInnes { 49537820f0SBarry Smith int bs = 1; /* Default block size = 1 */ 50537820f0SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 513b2fbd54SBarry Smith ierr = MatCreateMPIBDiag(mat->comm,PETSC_DECIDE,m,n,0,bs,PETSC_NULL, 52b4fd4287SBarry Smith PETSC_NULL,M); CHKERRQ(ierr); 538b6375c0SLois Curfman McInnes break; 548b6375c0SLois Curfman McInnes } 5594a9d846SBarry Smith case MATSEQBAIJ: 5694a9d846SBarry Smith ierr = MatCreateSeqBAIJ(mat->comm,1,m,n,0,PETSC_NULL,M); CHKERRQ(ierr); 5794a9d846SBarry Smith break; 5894a9d846SBarry Smith case MATMPIBAIJ: 5994a9d846SBarry Smith ierr = MatCreateMPIBAIJ(mat->comm,1,PETSC_DECIDE,PETSC_DECIDE, 6094a9d846SBarry Smith m,n,0,PETSC_NULL,0,PETSC_NULL,M); CHKERRQ(ierr); 6194a9d846SBarry Smith break; 628b6375c0SLois Curfman McInnes default: 63e3372554SBarry Smith SETERRQ(1,0,"Matrix type is not currently supported"); 648b6375c0SLois Curfman McInnes } 658b6375c0SLois Curfman McInnes ierr = MatGetOwnershipRange(*M,&rstart,&rend); CHKERRQ(ierr); 668b6375c0SLois Curfman McInnes for (i=rstart; i<rend; i++) { 678b6375c0SLois Curfman McInnes ierr = MatGetRow(mat,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 688b6375c0SLois Curfman McInnes ierr = MatSetValues(*M,1,&i,nz,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 698b6375c0SLois Curfman McInnes ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 708b6375c0SLois Curfman McInnes } 716d4a8577SBarry Smith ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 726d4a8577SBarry Smith ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 738b6375c0SLois Curfman McInnes return 0; 748b6375c0SLois Curfman McInnes } 75