xref: /petsc/src/mat/tests/ex218.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests MatShellTestMult()\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown typedef struct _n_User *User;
7*c4762a1bSJed Brown struct _n_User {
8*c4762a1bSJed Brown   Mat B;
9*c4762a1bSJed Brown };
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown static PetscErrorCode MatMult_User(Mat A,Vec X,Vec Y)
12*c4762a1bSJed Brown {
13*c4762a1bSJed Brown   User           user;
14*c4762a1bSJed Brown   PetscErrorCode ierr;
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown   PetscFunctionBegin;
17*c4762a1bSJed Brown   ierr = MatShellGetContext(A,&user);CHKERRQ(ierr);
18*c4762a1bSJed Brown   ierr = MatMult(user->B,X,Y);CHKERRQ(ierr);
19*c4762a1bSJed Brown   PetscFunctionReturn(0);
20*c4762a1bSJed Brown }
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown static PetscErrorCode MatMultTranspose_User(Mat A,Vec X,Vec Y)
23*c4762a1bSJed Brown {
24*c4762a1bSJed Brown   User           user;
25*c4762a1bSJed Brown   PetscErrorCode ierr;
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown   PetscFunctionBegin;
28*c4762a1bSJed Brown   ierr = MatShellGetContext(A,&user);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = MatMultTranspose(user->B,X,Y);CHKERRQ(ierr);
30*c4762a1bSJed Brown   PetscFunctionReturn(0);
31*c4762a1bSJed Brown }
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown static PetscErrorCode MyFunction(void *ctx,Vec x,Vec y)
34*c4762a1bSJed Brown {
35*c4762a1bSJed Brown   User           user = (User) ctx;
36*c4762a1bSJed Brown   PetscErrorCode ierr;
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown   PetscFunctionBegin;
39*c4762a1bSJed Brown   ierr = MatMult(user->B,x,y);CHKERRQ(ierr);
40*c4762a1bSJed Brown   PetscFunctionReturn(0);
41*c4762a1bSJed Brown }
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown int main(int argc,char **args)
44*c4762a1bSJed Brown {
45*c4762a1bSJed Brown   const PetscInt    inds[]  = {0,1};
46*c4762a1bSJed Brown   PetscScalar       avals[] = {2,3,5,7};
47*c4762a1bSJed Brown   Mat               S;
48*c4762a1bSJed Brown   User              user;
49*c4762a1bSJed Brown   PetscErrorCode    ierr;
50*c4762a1bSJed Brown   Vec               base;
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
53*c4762a1bSJed Brown   ierr = PetscNew(&user);CHKERRQ(ierr);
54*c4762a1bSJed Brown   ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD,2,2,2,NULL,&user->B);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = MatSetUp(user->B);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = MatSetValues(user->B,2,inds,2,inds,avals,INSERT_VALUES);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = MatAssemblyBegin(user->B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr = MatAssemblyEnd(user->B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = MatCreateVecs(user->B,&base,NULL);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = MatCreateShell(PETSC_COMM_WORLD,2,2,2,2,user,&S);CHKERRQ(ierr);
61*c4762a1bSJed Brown   ierr = MatSetUp(S);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MatMult_User);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = MatShellSetOperation(S,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_User);CHKERRQ(ierr);
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown   ierr = MatShellTestMult(S,MyFunction,base,user,NULL);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = MatShellTestMultTranspose(S,MyFunction,base,user,NULL);CHKERRQ(ierr);
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown   ierr = VecDestroy(&base);CHKERRQ(ierr);
69*c4762a1bSJed Brown   ierr = MatDestroy(&user->B);CHKERRQ(ierr);
70*c4762a1bSJed Brown   ierr = MatDestroy(&S);CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = PetscFree(user);CHKERRQ(ierr);
72*c4762a1bSJed Brown   ierr = PetscFinalize();
73*c4762a1bSJed Brown   return ierr;
74*c4762a1bSJed Brown }
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown /*TEST
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown    test:
80*c4762a1bSJed Brown      args: -mat_shell_test_mult_view -mat_shell_test_mult_transpose_view
81*c4762a1bSJed Brown 
82*c4762a1bSJed Brown TEST*/
83