xref: /petsc/src/mat/tests/ex175.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatCreateHermitianTranspose().\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7*d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat         C, C_htransposed, Cht, C_empty;
9c4762a1bSJed Brown   PetscInt    i, j, m = 10, n = 10;
10c4762a1bSJed Brown   PetscScalar v;
11c4762a1bSJed Brown 
12327415f7SBarry Smith   PetscFunctionBeginUser;
139566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
14c4762a1bSJed Brown   /* Create a complex non-hermitian matrix */
159566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &C));
169566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m, n));
179566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
189566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
19c4762a1bSJed Brown   for (i = 0; i < m; i++) {
20c4762a1bSJed Brown     for (j = 0; j < n; j++) {
21c4762a1bSJed Brown       v = 0.0 - 1.0 * PETSC_i;
229566063dSJacob Faibussowitsch       if (i > j && i - j < 2) PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
23c4762a1bSJed Brown     }
24c4762a1bSJed Brown   }
259566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
27c4762a1bSJed Brown 
289566063dSJacob Faibussowitsch   PetscCall(MatCreateHermitianTranspose(C, &C_htransposed));
29c4762a1bSJed Brown 
309566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_SELF));
319566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(C_htransposed, MAT_COPY_VALUES, &Cht));
329566063dSJacob Faibussowitsch   PetscCall(MatView(Cht, PETSC_VIEWER_STDOUT_SELF));
339566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(C_htransposed, MAT_DO_NOT_COPY_VALUES, &C_empty));
349566063dSJacob Faibussowitsch   PetscCall(MatView(C_empty, PETSC_VIEWER_STDOUT_SELF));
35c4762a1bSJed Brown 
369566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
379566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C_htransposed));
389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Cht));
399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C_empty));
409566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
41b122ec5aSJacob Faibussowitsch   return 0;
42c4762a1bSJed Brown }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown /*TEST
45c4762a1bSJed Brown 
46c4762a1bSJed Brown    build:
47c4762a1bSJed Brown      requires: complex
48c4762a1bSJed Brown 
49c4762a1bSJed Brown    test:
50c4762a1bSJed Brown      output_file: output/ex175.out
51c4762a1bSJed Brown 
52c4762a1bSJed Brown TEST*/
53