1af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 26098dc10SBarry Smith 3b3d09e86SJed Brown PetscErrorCode MatConvert_Shell(Mat oldmat, MatType newtype,MatReuse reuse,Mat *newmat) 4dfbe8321SBarry Smith { 5676c34cdSKris Buschelman Mat mat; 66098dc10SBarry Smith Vec in,out; 7dfbe8321SBarry Smith PetscErrorCode ierr; 86696f472SJed Brown PetscInt i,m,n,M,N,*rows,start; 96098dc10SBarry Smith MPI_Comm comm; 106696f472SJed Brown PetscScalar *array; 116098dc10SBarry Smith 123a40ed3dSBarry Smith PetscFunctionBegin; 13ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)oldmat,&comm);CHKERRQ(ierr); 146098dc10SBarry Smith 156696f472SJed Brown ierr = MatGetOwnershipRange(oldmat,&start,NULL);CHKERRQ(ierr); 166696f472SJed Brown ierr = MatCreateVecs(oldmat,&in,&out);CHKERRQ(ierr); 176696f472SJed Brown ierr = MatGetLocalSize(oldmat,&m,&n);CHKERRQ(ierr); 186696f472SJed Brown ierr = MatGetSize(oldmat,&M,&N);CHKERRQ(ierr); 196696f472SJed Brown ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 202205254eSKarl Rupp for (i=0; i<m; i++) rows[i] = start + i; 216098dc10SBarry Smith 22f69a0ea3SMatthew Knepley ierr = MatCreate(comm,&mat);CHKERRQ(ierr); 236696f472SJed Brown ierr = MatSetSizes(mat,m,n,M,N);CHKERRQ(ierr); 244a6a16e5SBarry Smith ierr = MatSetType(mat,newtype);CHKERRQ(ierr); 2533d57670SJed Brown ierr = MatSetBlockSizesFromMats(mat,oldmat,oldmat);CHKERRQ(ierr); 264a6a16e5SBarry Smith ierr = MatSetUp(mat);CHKERRQ(ierr); 276696f472SJed Brown ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 286098dc10SBarry Smith 296696f472SJed Brown for (i=0; i<N; i++) { 306696f472SJed Brown ierr = VecZeroEntries(in);CHKERRQ(ierr); 316696f472SJed Brown ierr = VecSetValue(in,i,1.,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); 436bf464f9SBarry Smith ierr = VecDestroy(&in);CHKERRQ(ierr); 446bf464f9SBarry 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 48511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 4928be2f97SBarry Smith ierr = MatHeaderReplace(oldmat,&mat);CHKERRQ(ierr); 50c3d102feSKris Buschelman } else { 51676c34cdSKris Buschelman *newmat = mat; 52c3d102feSKris Buschelman } 533a40ed3dSBarry Smith PetscFunctionReturn(0); 546098dc10SBarry Smith } 55*251fad3fSStefano Zampini 56*251fad3fSStefano Zampini static PetscErrorCode MatGetDiagonal_CF(Mat A,Vec X) 57*251fad3fSStefano Zampini { 58*251fad3fSStefano Zampini Mat B; 59*251fad3fSStefano Zampini PetscErrorCode ierr; 60*251fad3fSStefano Zampini 61*251fad3fSStefano Zampini PetscFunctionBegin; 62*251fad3fSStefano Zampini ierr = MatShellGetContext(A,&B);CHKERRQ(ierr); 63*251fad3fSStefano Zampini if (!B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing user matrix"); 64*251fad3fSStefano Zampini ierr = MatGetDiagonal(B,X);CHKERRQ(ierr); 65*251fad3fSStefano Zampini PetscFunctionReturn(0); 66*251fad3fSStefano Zampini } 67*251fad3fSStefano Zampini 68*251fad3fSStefano Zampini static PetscErrorCode MatMult_CF(Mat A,Vec X,Vec Y) 69*251fad3fSStefano Zampini { 70*251fad3fSStefano Zampini Mat B; 71*251fad3fSStefano Zampini PetscErrorCode ierr; 72*251fad3fSStefano Zampini 73*251fad3fSStefano Zampini PetscFunctionBegin; 74*251fad3fSStefano Zampini ierr = MatShellGetContext(A,&B);CHKERRQ(ierr); 75*251fad3fSStefano Zampini if (!B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing user matrix"); 76*251fad3fSStefano Zampini ierr = MatMult(B,X,Y);CHKERRQ(ierr); 77*251fad3fSStefano Zampini PetscFunctionReturn(0); 78*251fad3fSStefano Zampini } 79*251fad3fSStefano Zampini 80*251fad3fSStefano Zampini static PetscErrorCode MatMultTranspose_CF(Mat A,Vec X,Vec Y) 81*251fad3fSStefano Zampini { 82*251fad3fSStefano Zampini Mat B; 83*251fad3fSStefano Zampini PetscErrorCode ierr; 84*251fad3fSStefano Zampini 85*251fad3fSStefano Zampini PetscFunctionBegin; 86*251fad3fSStefano Zampini ierr = MatShellGetContext(A,&B);CHKERRQ(ierr); 87*251fad3fSStefano Zampini if (!B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing user matrix"); 88*251fad3fSStefano Zampini ierr = MatMultTranspose(B,X,Y);CHKERRQ(ierr); 89*251fad3fSStefano Zampini PetscFunctionReturn(0); 90*251fad3fSStefano Zampini } 91*251fad3fSStefano Zampini 92*251fad3fSStefano Zampini static PetscErrorCode MatDestroy_CF(Mat A) 93*251fad3fSStefano Zampini { 94*251fad3fSStefano Zampini Mat B; 95*251fad3fSStefano Zampini PetscErrorCode ierr; 96*251fad3fSStefano Zampini 97*251fad3fSStefano Zampini PetscFunctionBegin; 98*251fad3fSStefano Zampini ierr = MatShellGetContext(A,&B);CHKERRQ(ierr); 99*251fad3fSStefano Zampini if (!B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing user matrix"); 100*251fad3fSStefano Zampini ierr = MatDestroy(&B);CHKERRQ(ierr); 101*251fad3fSStefano Zampini PetscFunctionReturn(0); 102*251fad3fSStefano Zampini } 103*251fad3fSStefano Zampini 104*251fad3fSStefano Zampini PetscErrorCode MatConvertFrom_Shell(Mat A, MatType newtype,MatReuse reuse,Mat *B) 105*251fad3fSStefano Zampini { 106*251fad3fSStefano Zampini Mat M; 107*251fad3fSStefano Zampini PetscBool flg; 108*251fad3fSStefano Zampini PetscErrorCode ierr; 109*251fad3fSStefano Zampini 110*251fad3fSStefano Zampini PetscFunctionBegin; 111*251fad3fSStefano Zampini ierr = PetscStrcmp(newtype,MATSHELL,&flg);CHKERRQ(ierr); 112*251fad3fSStefano Zampini if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only conversion to MATSHELL"); 113*251fad3fSStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 114*251fad3fSStefano Zampini ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 115*251fad3fSStefano Zampini ierr = MatCreateShell(PetscObjectComm((PetscObject)A),A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N,A,&M);CHKERRQ(ierr); 116*251fad3fSStefano Zampini ierr = MatShellSetOperation(M,MATOP_MULT, (void (*)(void))MatMult_CF);CHKERRQ(ierr); 117*251fad3fSStefano Zampini ierr = MatShellSetOperation(M,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_CF);CHKERRQ(ierr); 118*251fad3fSStefano Zampini ierr = MatShellSetOperation(M,MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_CF);CHKERRQ(ierr); 119*251fad3fSStefano Zampini ierr = MatShellSetOperation(M,MATOP_DESTROY, (void (*)(void))MatDestroy_CF);CHKERRQ(ierr); 120*251fad3fSStefano Zampini *B = M; 121*251fad3fSStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented"); 122*251fad3fSStefano Zampini PetscFunctionReturn(0); 123*251fad3fSStefano Zampini } 124