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