xref: /petsc/src/mat/tests/ex218.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatShellTestMult()\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown typedef struct _n_User *User;
7c4762a1bSJed Brown struct _n_User {
8c4762a1bSJed Brown   Mat B;
9c4762a1bSJed Brown };
10c4762a1bSJed Brown 
11c4762a1bSJed Brown static PetscErrorCode MatMult_User(Mat A,Vec X,Vec Y)
12c4762a1bSJed Brown {
13c4762a1bSJed Brown   User           user;
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   PetscFunctionBegin;
165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(A,&user));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->B,X,Y));
18c4762a1bSJed Brown   PetscFunctionReturn(0);
19c4762a1bSJed Brown }
20c4762a1bSJed Brown 
21c4762a1bSJed Brown static PetscErrorCode MatMultTranspose_User(Mat A,Vec X,Vec Y)
22c4762a1bSJed Brown {
23c4762a1bSJed Brown   User           user;
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   PetscFunctionBegin;
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(A,&user));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(user->B,X,Y));
28c4762a1bSJed Brown   PetscFunctionReturn(0);
29c4762a1bSJed Brown }
30c4762a1bSJed Brown 
31c4762a1bSJed Brown static PetscErrorCode MyFunction(void *ctx,Vec x,Vec y)
32c4762a1bSJed Brown {
33c4762a1bSJed Brown   User           user = (User) ctx;
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   PetscFunctionBegin;
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->B,x,y));
37c4762a1bSJed Brown   PetscFunctionReturn(0);
38c4762a1bSJed Brown }
39c4762a1bSJed Brown 
40c4762a1bSJed Brown int main(int argc,char **args)
41c4762a1bSJed Brown {
42c4762a1bSJed Brown   const PetscInt    inds[]  = {0,1};
43c4762a1bSJed Brown   PetscScalar       avals[] = {2,3,5,7};
44c4762a1bSJed Brown   Mat               S;
45c4762a1bSJed Brown   User              user;
46c4762a1bSJed Brown   Vec               base;
47c4762a1bSJed Brown 
48*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&user));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD,2,2,2,NULL,&user->B));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user->B));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(user->B,2,inds,2,inds,avals,INSERT_VALUES));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(user->B,MAT_FINAL_ASSEMBLY));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(user->B,MAT_FINAL_ASSEMBLY));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user->B,&base,NULL));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,2,2,2,2,user,&S));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(S));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MatMult_User));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(S,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_User));
60c4762a1bSJed Brown 
615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellTestMult(S,MyFunction,base,user,NULL));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellTestMultTranspose(S,MyFunction,base,user,NULL));
63c4762a1bSJed Brown 
645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&base));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->B));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&S));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user));
68*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
69*b122ec5aSJacob Faibussowitsch   return 0;
70c4762a1bSJed Brown }
71c4762a1bSJed Brown 
72c4762a1bSJed Brown /*TEST
73c4762a1bSJed Brown 
74c4762a1bSJed Brown    test:
75c4762a1bSJed Brown      args: -mat_shell_test_mult_view -mat_shell_test_mult_transpose_view
76c4762a1bSJed Brown 
77c4762a1bSJed Brown TEST*/
78