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