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; 76098dc10SBarry Smith MPI_Comm comm; 86696f472SJed Brown PetscScalar *array; 9*41319c1dSStefano Zampini PetscInt *dnnz,*onnz,*dnnzu,*onnzu; 10*41319c1dSStefano Zampini PetscInt cst,Nbs,mbs,nbs,rbs,cbs; 11*41319c1dSStefano Zampini PetscInt im,i,m,n,M,N,*rows,start; 12*41319c1dSStefano Zampini PetscErrorCode ierr; 136098dc10SBarry Smith 143a40ed3dSBarry Smith PetscFunctionBegin; 15ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)oldmat,&comm);CHKERRQ(ierr); 166098dc10SBarry Smith 176696f472SJed Brown ierr = MatGetOwnershipRange(oldmat,&start,NULL);CHKERRQ(ierr); 18*41319c1dSStefano Zampini ierr = MatGetOwnershipRangeColumn(oldmat,&cst,NULL);CHKERRQ(ierr); 196696f472SJed Brown ierr = MatCreateVecs(oldmat,&in,&out);CHKERRQ(ierr); 206696f472SJed Brown ierr = MatGetLocalSize(oldmat,&m,&n);CHKERRQ(ierr); 216696f472SJed Brown ierr = MatGetSize(oldmat,&M,&N);CHKERRQ(ierr); 226696f472SJed Brown ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 236098dc10SBarry Smith 24f69a0ea3SMatthew Knepley ierr = MatCreate(comm,&mat);CHKERRQ(ierr); 256696f472SJed Brown ierr = MatSetSizes(mat,m,n,M,N);CHKERRQ(ierr); 264a6a16e5SBarry Smith ierr = MatSetType(mat,newtype);CHKERRQ(ierr); 2733d57670SJed Brown ierr = MatSetBlockSizesFromMats(mat,oldmat,oldmat);CHKERRQ(ierr); 28*41319c1dSStefano Zampini ierr = MatGetBlockSizes(mat,&rbs,&cbs);CHKERRQ(ierr); 29*41319c1dSStefano Zampini mbs = m/rbs; 30*41319c1dSStefano Zampini nbs = n/cbs; 31*41319c1dSStefano Zampini Nbs = N/cbs; 32*41319c1dSStefano Zampini cst = cst/cbs; 33*41319c1dSStefano Zampini ierr = PetscMalloc4(mbs,&dnnz,mbs,&onnz,mbs,&dnnzu,mbs,&onnzu);CHKERRQ(ierr); 34*41319c1dSStefano Zampini for (i=0; i<mbs; i++) { 35*41319c1dSStefano Zampini dnnz[i] = nbs; 36*41319c1dSStefano Zampini onnz[i] = Nbs - nbs; 37*41319c1dSStefano Zampini dnnzu[i] = PetscMax(nbs - i,0); 38*41319c1dSStefano Zampini onnzu[i] = PetscMax(Nbs - (cst + nbs),0); 39*41319c1dSStefano Zampini } 40*41319c1dSStefano Zampini ierr = MatXAIJSetPreallocation(mat,PETSC_DECIDE,dnnz,onnz,dnnzu,onnzu);CHKERRQ(ierr); 41*41319c1dSStefano Zampini ierr = PetscFree4(dnnz,onnz,dnnzu,onnzu);CHKERRQ(ierr); 426696f472SJed Brown ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 43*41319c1dSStefano Zampini ierr = MatSetUp(mat);CHKERRQ(ierr); 446696f472SJed Brown for (i=0; i<N; i++) { 45*41319c1dSStefano Zampini PetscInt j; 46*41319c1dSStefano Zampini 476696f472SJed Brown ierr = VecZeroEntries(in);CHKERRQ(ierr); 486696f472SJed Brown ierr = VecSetValue(in,i,1.,INSERT_VALUES);CHKERRQ(ierr); 496098dc10SBarry Smith ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 506098dc10SBarry Smith ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 516098dc10SBarry Smith ierr = MatMult(oldmat,in,out);CHKERRQ(ierr); 526098dc10SBarry Smith ierr = VecGetArray(out,&array);CHKERRQ(ierr); 53*41319c1dSStefano Zampini for (j=0, im = 0; j<m; j++) { 54*41319c1dSStefano Zampini if (PetscAbsScalar(array[j]) == 0.0) continue; 55*41319c1dSStefano Zampini rows[im] = j+start; 56*41319c1dSStefano Zampini array[im] = array[j]; 57*41319c1dSStefano Zampini im++; 58*41319c1dSStefano Zampini } 59*41319c1dSStefano Zampini ierr = MatSetValues(mat,im,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 606098dc10SBarry Smith ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 616098dc10SBarry Smith } 62606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 636bf464f9SBarry Smith ierr = VecDestroy(&in);CHKERRQ(ierr); 646bf464f9SBarry Smith ierr = VecDestroy(&out);CHKERRQ(ierr); 65676c34cdSKris Buschelman ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 66676c34cdSKris Buschelman ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 67511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 6828be2f97SBarry Smith ierr = MatHeaderReplace(oldmat,&mat);CHKERRQ(ierr); 69c3d102feSKris Buschelman } else { 70676c34cdSKris Buschelman *newmat = mat; 71c3d102feSKris Buschelman } 723a40ed3dSBarry Smith PetscFunctionReturn(0); 736098dc10SBarry Smith } 74251fad3fSStefano Zampini 75251fad3fSStefano Zampini static PetscErrorCode MatGetDiagonal_CF(Mat A,Vec X) 76251fad3fSStefano Zampini { 77251fad3fSStefano Zampini Mat B; 78251fad3fSStefano Zampini PetscErrorCode ierr; 79251fad3fSStefano Zampini 80251fad3fSStefano Zampini PetscFunctionBegin; 81251fad3fSStefano Zampini ierr = MatShellGetContext(A,&B);CHKERRQ(ierr); 82251fad3fSStefano Zampini if (!B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing user matrix"); 83251fad3fSStefano Zampini ierr = MatGetDiagonal(B,X);CHKERRQ(ierr); 84251fad3fSStefano Zampini PetscFunctionReturn(0); 85251fad3fSStefano Zampini } 86251fad3fSStefano Zampini 87251fad3fSStefano Zampini static PetscErrorCode MatMult_CF(Mat A,Vec X,Vec Y) 88251fad3fSStefano Zampini { 89251fad3fSStefano Zampini Mat B; 90251fad3fSStefano Zampini PetscErrorCode ierr; 91251fad3fSStefano Zampini 92251fad3fSStefano Zampini PetscFunctionBegin; 93251fad3fSStefano Zampini ierr = MatShellGetContext(A,&B);CHKERRQ(ierr); 94251fad3fSStefano Zampini if (!B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing user matrix"); 95251fad3fSStefano Zampini ierr = MatMult(B,X,Y);CHKERRQ(ierr); 96251fad3fSStefano Zampini PetscFunctionReturn(0); 97251fad3fSStefano Zampini } 98251fad3fSStefano Zampini 99251fad3fSStefano Zampini static PetscErrorCode MatMultTranspose_CF(Mat A,Vec X,Vec Y) 100251fad3fSStefano Zampini { 101251fad3fSStefano Zampini Mat B; 102251fad3fSStefano Zampini PetscErrorCode ierr; 103251fad3fSStefano Zampini 104251fad3fSStefano Zampini PetscFunctionBegin; 105251fad3fSStefano Zampini ierr = MatShellGetContext(A,&B);CHKERRQ(ierr); 106251fad3fSStefano Zampini if (!B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing user matrix"); 107251fad3fSStefano Zampini ierr = MatMultTranspose(B,X,Y);CHKERRQ(ierr); 108251fad3fSStefano Zampini PetscFunctionReturn(0); 109251fad3fSStefano Zampini } 110251fad3fSStefano Zampini 111251fad3fSStefano Zampini static PetscErrorCode MatDestroy_CF(Mat A) 112251fad3fSStefano Zampini { 113251fad3fSStefano Zampini Mat B; 114251fad3fSStefano Zampini PetscErrorCode ierr; 115251fad3fSStefano Zampini 116251fad3fSStefano Zampini PetscFunctionBegin; 117251fad3fSStefano Zampini ierr = MatShellGetContext(A,&B);CHKERRQ(ierr); 118251fad3fSStefano Zampini if (!B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing user matrix"); 119251fad3fSStefano Zampini ierr = MatDestroy(&B);CHKERRQ(ierr); 120251fad3fSStefano Zampini PetscFunctionReturn(0); 121251fad3fSStefano Zampini } 122251fad3fSStefano Zampini 123251fad3fSStefano Zampini PetscErrorCode MatConvertFrom_Shell(Mat A, MatType newtype,MatReuse reuse,Mat *B) 124251fad3fSStefano Zampini { 125251fad3fSStefano Zampini Mat M; 126251fad3fSStefano Zampini PetscBool flg; 127251fad3fSStefano Zampini PetscErrorCode ierr; 128251fad3fSStefano Zampini 129251fad3fSStefano Zampini PetscFunctionBegin; 130251fad3fSStefano Zampini ierr = PetscStrcmp(newtype,MATSHELL,&flg);CHKERRQ(ierr); 131251fad3fSStefano Zampini if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only conversion to MATSHELL"); 132251fad3fSStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 133251fad3fSStefano Zampini ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 134251fad3fSStefano Zampini ierr = MatCreateShell(PetscObjectComm((PetscObject)A),A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N,A,&M);CHKERRQ(ierr); 1351d9d9a34SStefano Zampini ierr = MatSetBlockSizesFromMats(M,A,A);CHKERRQ(ierr); 136251fad3fSStefano Zampini ierr = MatShellSetOperation(M,MATOP_MULT, (void (*)(void))MatMult_CF);CHKERRQ(ierr); 137251fad3fSStefano Zampini ierr = MatShellSetOperation(M,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_CF);CHKERRQ(ierr); 138251fad3fSStefano Zampini ierr = MatShellSetOperation(M,MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_CF);CHKERRQ(ierr); 139251fad3fSStefano Zampini ierr = MatShellSetOperation(M,MATOP_DESTROY, (void (*)(void))MatDestroy_CF);CHKERRQ(ierr); 140251fad3fSStefano Zampini *B = M; 141251fad3fSStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented"); 142251fad3fSStefano Zampini PetscFunctionReturn(0); 143251fad3fSStefano Zampini } 144