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