xref: /petsc/src/mat/tests/ex235.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
105f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(S,&A));
115f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(S,&A));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&ll,&rr));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&d,NULL));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(ll,NULL));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(rr,NULL));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(d,NULL));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(A,3.0));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(A,-4.0));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(A,8.0));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDiagonalSet(A,d,ADD_VALUES));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(A,9.0));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(A,8.0));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(ll,NULL));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(rr,NULL));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDiagonalScale(A,ll,rr));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(A,2.0));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(A,11.0));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(d,NULL));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDiagonalSet(A,d,ADD_VALUES));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(ll,NULL));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(rr,NULL));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDiagonalScale(A,ll,rr));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(A,5.0));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(A,7.0));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetDiagonal(A,d));
54c4762a1bSJed Brown   *D   = d;
555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ll));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
67*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
69c4762a1bSJed Brown 
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,m,7,NULL,6,NULL,&Aij));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(Aij,NULL));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(Aij,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE));
73c4762a1bSJed Brown 
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,m,Aij,&A));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(A,MATOP_MULT,(void (*)(void)) myMult));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void (*)(void)) myGetDiagonal));
77c4762a1bSJed Brown 
785f80ce2aSJacob Faibussowitsch   CHKERRQ(shiftandscale(A,&Adiag));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatComputeOperator(A,NULL,&B));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(shiftandscale(Aij,&Aijdiag));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(Aij,-1.0,B,DIFFERENT_NONZERO_PATTERN));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(Aij,NORM_FROBENIUS,&Aijnorm));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(B,NORM_FROBENIUS,&Bnorm));
842c71b3e2SJacob 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));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(Aijdiag,-1.0,Adiag));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(Adiag,NORM_2,&dnorm));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(Aijdiag,NORM_2,&Aijdiagnorm));
882c71b3e2SJacob 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));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Aij));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Adiag));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Aijdiag));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
94*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
95*b122ec5aSJacob 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