1*c4762a1bSJed Brown static char help[] = "Test combinatations of scalings, shifts and get diagonal of MATSHELL\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscmat.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown static PetscErrorCode myMult(Mat S,Vec x,Vec y) 6*c4762a1bSJed Brown { 7*c4762a1bSJed Brown Mat A; 8*c4762a1bSJed Brown PetscErrorCode ierr; 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown PetscFunctionBegin; 11*c4762a1bSJed Brown ierr = MatShellGetContext(S,&A);CHKERRQ(ierr); 12*c4762a1bSJed Brown ierr = MatMult(A,x,y);CHKERRQ(ierr); 13*c4762a1bSJed Brown PetscFunctionReturn(0); 14*c4762a1bSJed Brown } 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown static PetscErrorCode myGetDiagonal(Mat S,Vec d) 17*c4762a1bSJed Brown { 18*c4762a1bSJed Brown Mat A; 19*c4762a1bSJed Brown PetscErrorCode ierr; 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown PetscFunctionBegin; 22*c4762a1bSJed Brown ierr = MatShellGetContext(S,&A);CHKERRQ(ierr); 23*c4762a1bSJed Brown ierr = MatGetDiagonal(A,d);CHKERRQ(ierr); 24*c4762a1bSJed Brown PetscFunctionReturn(0); 25*c4762a1bSJed Brown } 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown static PetscErrorCode shiftandscale(Mat A,Vec *D) 28*c4762a1bSJed Brown { 29*c4762a1bSJed Brown PetscErrorCode ierr; 30*c4762a1bSJed Brown Vec ll,d,rr; 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown PetscFunctionBegin; 33*c4762a1bSJed Brown ierr = MatCreateVecs(A,&ll,&rr);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = MatCreateVecs(A,&d,NULL);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = VecSetRandom(ll,NULL);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = VecSetRandom(rr,NULL);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = VecSetRandom(d,NULL);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = MatScale(A,3.0);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = MatShift(A,-4.0);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = MatScale(A,8.0);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = MatDiagonalSet(A,d,ADD_VALUES);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = MatShift(A,9.0);CHKERRQ(ierr); 43*c4762a1bSJed Brown ierr = MatScale(A,8.0);CHKERRQ(ierr); 44*c4762a1bSJed Brown ierr = VecSetRandom(ll,NULL);CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = VecSetRandom(rr,NULL);CHKERRQ(ierr); 46*c4762a1bSJed Brown ierr = MatDiagonalScale(A,ll,rr);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = MatShift(A,2.0);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = MatScale(A,11.0);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = VecSetRandom(d,NULL);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = MatDiagonalSet(A,d,ADD_VALUES);CHKERRQ(ierr); 51*c4762a1bSJed Brown ierr = VecSetRandom(ll,NULL);CHKERRQ(ierr); 52*c4762a1bSJed Brown ierr = VecSetRandom(rr,NULL);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = MatDiagonalScale(A,ll,rr);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = MatShift(A,5.0);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = MatScale(A,7.0);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = MatGetDiagonal(A,d);CHKERRQ(ierr); 57*c4762a1bSJed Brown *D = d; 58*c4762a1bSJed Brown ierr = VecDestroy(&ll);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = VecDestroy(&rr);CHKERRQ(ierr); 60*c4762a1bSJed Brown PetscFunctionReturn(0); 61*c4762a1bSJed Brown } 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown int main(int argc,char **args) 64*c4762a1bSJed Brown { 65*c4762a1bSJed Brown Mat A,Aij,B; 66*c4762a1bSJed Brown Vec Adiag,Aijdiag; 67*c4762a1bSJed Brown PetscErrorCode ierr; 68*c4762a1bSJed Brown PetscInt m = 3; 69*c4762a1bSJed Brown PetscReal Aijnorm,Aijdiagnorm,Bnorm,dnorm; 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 72*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,m,7,NULL,6,NULL,&Aij);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = MatSetRandom(Aij,NULL);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = MatSetOption(Aij,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown ierr = MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,m,Aij,&A);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = MatShellSetOperation(A,MATOP_MULT,(void (*)(void)) myMult);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void (*)(void)) myGetDiagonal);CHKERRQ(ierr); 81*c4762a1bSJed Brown 82*c4762a1bSJed Brown ierr = shiftandscale(A,&Adiag);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = MatComputeOperator(A,NULL,&B);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = shiftandscale(Aij,&Aijdiag);CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = MatAXPY(Aij,-1.0,B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = MatNorm(Aij,NORM_FROBENIUS,&Aijnorm);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = MatNorm(B,NORM_FROBENIUS,&Bnorm);CHKERRQ(ierr); 88*c4762a1bSJed Brown if (Aijnorm/Bnorm > 100.0*PETSC_MACHINE_EPSILON) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Altered matrices do not match, norm of difference %g",Aijnorm/Bnorm); 89*c4762a1bSJed Brown ierr = VecAXPY(Aijdiag,-1.0,Adiag);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = VecNorm(Adiag,NORM_2,&dnorm);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = VecNorm(Aijdiag,NORM_2,&Aijdiagnorm);CHKERRQ(ierr); 92*c4762a1bSJed Brown if (Aijdiagnorm/dnorm > 100.0*PETSC_MACHINE_EPSILON) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Altered matrices diagonals do not match, norm of difference %g",Aijdiagnorm/dnorm); 93*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 94*c4762a1bSJed Brown ierr = MatDestroy(&Aij);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = VecDestroy(&Adiag);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = VecDestroy(&Aijdiag);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = PetscFinalize(); 99*c4762a1bSJed Brown return ierr; 100*c4762a1bSJed Brown } 101*c4762a1bSJed Brown 102*c4762a1bSJed Brown /*TEST 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown test: 105*c4762a1bSJed Brown nsize: {{1 2 3 4}} 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown TEST*/ 108