1 static char help[] = "Tests MatCreateHermitianTranspose().\n\n";
2
3 #include <petscmat.h>
4
main(int argc,char ** args)5 int main(int argc, char **args)
6 {
7 Mat C, C_htransposed, Cht, C_empty, Cht_empty;
8 PetscInt i, j, m = 10, n = 10;
9 PetscScalar v;
10
11 PetscFunctionBeginUser;
12 PetscCall(PetscInitialize(&argc, &args, NULL, help));
13 /* Create a complex non-hermitian matrix */
14 PetscCall(MatCreate(PETSC_COMM_SELF, &C));
15 PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m, n));
16 PetscCall(MatSetFromOptions(C));
17 PetscCall(MatSetUp(C));
18 for (i = 0; i < m; i++) {
19 for (j = 0; j < n; j++) {
20 v = 0.0 - 1.0 * PETSC_i;
21 if (i > j && i - j < 2) PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
22 }
23 }
24 PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
25 PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
26
27 PetscCall(MatCreateHermitianTranspose(C, &C_htransposed));
28
29 PetscCall(MatView(C, PETSC_VIEWER_STDOUT_SELF));
30 PetscCall(MatDuplicate(C_htransposed, MAT_COPY_VALUES, &Cht));
31 PetscCall(MatConvert(Cht, MATAIJ, MAT_INPLACE_MATRIX, &Cht));
32 PetscCall(MatView(Cht, PETSC_VIEWER_STDOUT_SELF));
33 PetscCall(MatDuplicate(C_htransposed, MAT_DO_NOT_COPY_VALUES, &C_empty));
34 PetscCall(MatHermitianTransposeGetMat(C_empty, &Cht_empty));
35 PetscCall(MatHermitianTranspose(Cht_empty, MAT_INPLACE_MATRIX, &Cht_empty));
36 PetscCall(MatView(Cht_empty, PETSC_VIEWER_STDOUT_SELF));
37
38 PetscCall(MatDestroy(&C));
39 PetscCall(MatDestroy(&C_htransposed));
40 PetscCall(MatDestroy(&Cht));
41 PetscCall(MatDestroy(&C_empty));
42 PetscCall(PetscFinalize());
43 return 0;
44 }
45
46 /*TEST
47
48 build:
49 requires: complex
50
51 test:
52 output_file: output/ex175.out
53
54 TEST*/
55