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