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