xref: /petsc/src/mat/tests/ex197.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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*327415f7SBarry Smith   PetscFunctionBeginUser;
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
159566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
169566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2));
179566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATAIJ));
189566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
199566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   i = 0; j = 0; v = 2.0;
229566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
23c4762a1bSJed Brown   i = 0; j = 1; v = 3.0 + 4.0*PETSC_i;
249566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
25c4762a1bSJed Brown   i = 1; j = 0; v = 5.0 + 6.0*PETSC_i;
269566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
27c4762a1bSJed Brown   i = 1; j = 1; v = 7.0 + 8.0*PETSC_i;
289566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /* Create vectors */
339566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&y));
349566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(y,PETSC_DECIDE,2));
359566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(y));
369566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(y,&ys));
379566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(y,&x));
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   i = 0; v = 10.0 + 11.0*PETSC_i;
409566063dSJacob Faibussowitsch   PetscCall(VecSetValues(x,1,&i,&v,INSERT_VALUES));
41c4762a1bSJed Brown   i = 1; v = 100.0 + 120.0*PETSC_i;
429566063dSJacob Faibussowitsch   PetscCall(VecSetValues(x,1,&i,&v,INSERT_VALUES));
439566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
449566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(A,x,y));
479566063dSJacob Faibussowitsch   PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
489566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTransposeAdd(A,x,y,ys));
499566063dSJacob Faibussowitsch   PetscCall(VecView(ys,PETSC_VIEWER_STDOUT_WORLD));
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&B));
529566063dSJacob Faibussowitsch   PetscCall(MatCreateHermitianTranspose(A,&C));
539566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTransposeEqual(B,C,4,&flg));
5428b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"B^Hx != C^Hx");
559566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTransposeAddEqual(B,C,4,&flg));
5628b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"y+B^Hx != y+C^Hx");
579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
589566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
59e573050aSPierre Jolivet 
609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ys));
659566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
66b122ec5aSJacob Faibussowitsch   return 0;
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