xref: /petsc/src/mat/tests/ex235.c (revision 2286efddd54511ab18e8e2adb1e023c4bf8f0b92)
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 
myMult(Mat S,Vec x,Vec y)5d71ae5a4SJacob Faibussowitsch static PetscErrorCode myMult(Mat S, Vec x, Vec y)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   Mat A;
8c4762a1bSJed Brown 
9c4762a1bSJed Brown   PetscFunctionBegin;
109566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(S, &A));
119566063dSJacob Faibussowitsch   PetscCall(MatMult(A, x, y));
123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13c4762a1bSJed Brown }
14c4762a1bSJed Brown 
myGetDiagonal(Mat S,Vec d)15d71ae5a4SJacob Faibussowitsch static PetscErrorCode myGetDiagonal(Mat S, Vec d)
16d71ae5a4SJacob Faibussowitsch {
17c4762a1bSJed Brown   Mat A;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(S, &A));
219566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(A, d));
223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23c4762a1bSJed Brown }
24c4762a1bSJed Brown 
shiftandscale(Mat A,Vec * D)25d71ae5a4SJacob Faibussowitsch static PetscErrorCode shiftandscale(Mat A, Vec *D)
26d71ae5a4SJacob Faibussowitsch {
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));
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58c4762a1bSJed Brown }
59c4762a1bSJed Brown 
main(int argc,char ** args)60d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
61d71ae5a4SJacob Faibussowitsch {
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 
67327415f7SBarry Smith   PetscFunctionBeginUser;
68c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, m, m, 7, NULL, 6, NULL, &Aij));
729566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(Aij, NULL));
739566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Aij, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
74c4762a1bSJed Brown 
759566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, m, m, Aij, &A));
76*57d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_MULT, (PetscErrorCodeFn *)myMult));
77*57d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)myGetDiagonal));
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch   PetscCall(shiftandscale(A, &Adiag));
809566063dSJacob Faibussowitsch   PetscCall(MatComputeOperator(A, NULL, &B));
819566063dSJacob Faibussowitsch   PetscCall(shiftandscale(Aij, &Aijdiag));
829566063dSJacob Faibussowitsch   PetscCall(MatAXPY(Aij, -1.0, B, DIFFERENT_NONZERO_PATTERN));
839566063dSJacob Faibussowitsch   PetscCall(MatNorm(Aij, NORM_FROBENIUS, &Aijnorm));
849566063dSJacob Faibussowitsch   PetscCall(MatNorm(B, NORM_FROBENIUS, &Bnorm));
8508401ef6SPierre 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));
869566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Aijdiag, -1.0, Adiag));
879566063dSJacob Faibussowitsch   PetscCall(VecNorm(Adiag, NORM_2, &dnorm));
889566063dSJacob Faibussowitsch   PetscCall(VecNorm(Aijdiag, NORM_2, &Aijdiagnorm));
8908401ef6SPierre 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));
909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Aij));
929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Adiag));
939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Aijdiag));
949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
959566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
96b122ec5aSJacob Faibussowitsch   return 0;
97c4762a1bSJed Brown }
98c4762a1bSJed Brown 
99c4762a1bSJed Brown /*TEST
100c4762a1bSJed Brown 
101c4762a1bSJed Brown     test:
102c4762a1bSJed Brown       nsize: {{1 2 3 4}}
1033886731fSPierre Jolivet       output_file: output/empty.out
104c4762a1bSJed Brown 
105c4762a1bSJed Brown TEST*/
106