xref: /petsc/src/mat/tests/ex197.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown static char help[] = "Test MatMultHermitianTranspose() and MatMultHermitianTransposeAdd().\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc,char **args)
6c4762a1bSJed Brown {
7e573050aSPierre Jolivet   Mat         A,B,C;
8c4762a1bSJed Brown   Vec         x,y,ys;
9c4762a1bSJed Brown   PetscInt    i,j;
10c4762a1bSJed Brown   PetscScalar v;
11e573050aSPierre Jolivet   PetscBool   flg;
12c4762a1bSJed Brown 
13*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
14*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
15*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2));
16*9566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATAIJ));
17*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
18*9566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   i = 0; j = 0; v = 2.0;
21*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
22c4762a1bSJed Brown   i = 0; j = 1; v = 3.0 + 4.0*PETSC_i;
23*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
24c4762a1bSJed Brown   i = 1; j = 0; v = 5.0 + 6.0*PETSC_i;
25*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
26c4762a1bSJed Brown   i = 1; j = 1; v = 7.0 + 8.0*PETSC_i;
27*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
28*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
29*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   /* Create vectors */
32*9566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&y));
33*9566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(y,PETSC_DECIDE,2));
34*9566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(y));
35*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(y,&ys));
36*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(y,&x));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   i = 0; v = 10.0 + 11.0*PETSC_i;
39*9566063dSJacob Faibussowitsch   PetscCall(VecSetValues(x,1,&i,&v,INSERT_VALUES));
40c4762a1bSJed Brown   i = 1; v = 100.0 + 120.0*PETSC_i;
41*9566063dSJacob Faibussowitsch   PetscCall(VecSetValues(x,1,&i,&v,INSERT_VALUES));
42*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
43*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
44c4762a1bSJed Brown 
45*9566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(A,x,y));
46*9566063dSJacob Faibussowitsch   PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
47*9566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTransposeAdd(A,x,y,ys));
48*9566063dSJacob Faibussowitsch   PetscCall(VecView(ys,PETSC_VIEWER_STDOUT_WORLD));
49c4762a1bSJed Brown 
50*9566063dSJacob Faibussowitsch   PetscCall(MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&B));
51*9566063dSJacob Faibussowitsch   PetscCall(MatCreateHermitianTranspose(A,&C));
52*9566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTransposeEqual(B,C,4,&flg));
5328b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"B^Hx != C^Hx");
54*9566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTransposeAddEqual(B,C,4,&flg));
5528b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"y+B^Hx != y+C^Hx");
56*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
57*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
58e573050aSPierre Jolivet 
59*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
60c4762a1bSJed Brown 
61*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
62*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
63*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ys));
64*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
65b122ec5aSJacob Faibussowitsch   return 0;
66c4762a1bSJed Brown }
67c4762a1bSJed Brown 
68c4762a1bSJed Brown /*TEST
69c4762a1bSJed Brown 
70c4762a1bSJed Brown    build:
71c4762a1bSJed Brown       requires: complex
72c4762a1bSJed Brown    test:
73c4762a1bSJed Brown 
74c4762a1bSJed Brown    test:
75c4762a1bSJed Brown       suffix: 2
76c4762a1bSJed Brown       nsize: 2
77c4762a1bSJed Brown 
78c4762a1bSJed Brown TEST*/
79