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 9c4762a1bSJed Brown PetscFunctionBegin; 109566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(S,&A)); 119566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,y)); 12c4762a1bSJed Brown PetscFunctionReturn(0); 13c4762a1bSJed Brown } 14c4762a1bSJed Brown 15c4762a1bSJed Brown static PetscErrorCode myGetDiagonal(Mat S,Vec d) 16c4762a1bSJed Brown { 17c4762a1bSJed Brown Mat A; 18c4762a1bSJed Brown 19c4762a1bSJed Brown PetscFunctionBegin; 209566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(S,&A)); 219566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A,d)); 22c4762a1bSJed Brown PetscFunctionReturn(0); 23c4762a1bSJed Brown } 24c4762a1bSJed Brown 25c4762a1bSJed Brown static PetscErrorCode shiftandscale(Mat A,Vec *D) 26c4762a1bSJed Brown { 27c4762a1bSJed Brown Vec ll,d,rr; 28c4762a1bSJed Brown 29c4762a1bSJed Brown PetscFunctionBegin; 309566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&ll,&rr)); 319566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&d,NULL)); 329566063dSJacob Faibussowitsch PetscCall(VecSetRandom(ll,NULL)); 339566063dSJacob Faibussowitsch PetscCall(VecSetRandom(rr,NULL)); 349566063dSJacob Faibussowitsch PetscCall(VecSetRandom(d,NULL)); 359566063dSJacob Faibussowitsch PetscCall(MatScale(A,3.0)); 369566063dSJacob Faibussowitsch PetscCall(MatShift(A,-4.0)); 379566063dSJacob Faibussowitsch PetscCall(MatScale(A,8.0)); 389566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A,d,ADD_VALUES)); 399566063dSJacob Faibussowitsch PetscCall(MatShift(A,9.0)); 409566063dSJacob Faibussowitsch PetscCall(MatScale(A,8.0)); 419566063dSJacob Faibussowitsch PetscCall(VecSetRandom(ll,NULL)); 429566063dSJacob Faibussowitsch PetscCall(VecSetRandom(rr,NULL)); 439566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A,ll,rr)); 449566063dSJacob Faibussowitsch PetscCall(MatShift(A,2.0)); 459566063dSJacob Faibussowitsch PetscCall(MatScale(A,11.0)); 469566063dSJacob Faibussowitsch PetscCall(VecSetRandom(d,NULL)); 479566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A,d,ADD_VALUES)); 489566063dSJacob Faibussowitsch PetscCall(VecSetRandom(ll,NULL)); 499566063dSJacob Faibussowitsch PetscCall(VecSetRandom(rr,NULL)); 509566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A,ll,rr)); 519566063dSJacob Faibussowitsch PetscCall(MatShift(A,5.0)); 529566063dSJacob Faibussowitsch PetscCall(MatScale(A,7.0)); 539566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A,d)); 54c4762a1bSJed Brown *D = d; 559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ll)); 569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rr)); 57c4762a1bSJed Brown PetscFunctionReturn(0); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown int main(int argc,char **args) 61c4762a1bSJed Brown { 62c4762a1bSJed Brown Mat A,Aij,B; 63c4762a1bSJed Brown Vec Adiag,Aijdiag; 64c4762a1bSJed Brown PetscInt m = 3; 65c4762a1bSJed Brown PetscReal Aijnorm,Aijdiagnorm,Bnorm,dnorm; 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 69c4762a1bSJed Brown 709566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,m,7,NULL,6,NULL,&Aij)); 719566063dSJacob Faibussowitsch PetscCall(MatSetRandom(Aij,NULL)); 729566063dSJacob Faibussowitsch PetscCall(MatSetOption(Aij,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE)); 73c4762a1bSJed Brown 749566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,m,Aij,&A)); 759566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A,MATOP_MULT,(void (*)(void)) myMult)); 769566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void (*)(void)) myGetDiagonal)); 77c4762a1bSJed Brown 789566063dSJacob Faibussowitsch PetscCall(shiftandscale(A,&Adiag)); 799566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(A,NULL,&B)); 809566063dSJacob Faibussowitsch PetscCall(shiftandscale(Aij,&Aijdiag)); 819566063dSJacob Faibussowitsch PetscCall(MatAXPY(Aij,-1.0,B,DIFFERENT_NONZERO_PATTERN)); 829566063dSJacob Faibussowitsch PetscCall(MatNorm(Aij,NORM_FROBENIUS,&Aijnorm)); 839566063dSJacob Faibussowitsch PetscCall(MatNorm(B,NORM_FROBENIUS,&Bnorm)); 84*08401ef6SPierre Jolivet PetscCheck(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)); 859566063dSJacob Faibussowitsch PetscCall(VecAXPY(Aijdiag,-1.0,Adiag)); 869566063dSJacob Faibussowitsch PetscCall(VecNorm(Adiag,NORM_2,&dnorm)); 879566063dSJacob Faibussowitsch PetscCall(VecNorm(Aijdiag,NORM_2,&Aijdiagnorm)); 88*08401ef6SPierre Jolivet PetscCheck(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)); 899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aij)); 919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Adiag)); 929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Aijdiag)); 939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 949566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 95b122ec5aSJacob Faibussowitsch return 0; 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98c4762a1bSJed Brown /*TEST 99c4762a1bSJed Brown 100c4762a1bSJed Brown test: 101c4762a1bSJed Brown nsize: {{1 2 3 4}} 102c4762a1bSJed Brown 103c4762a1bSJed Brown TEST*/ 104